149bd79ccSDebojyoti Ghosh 249bd79ccSDebojyoti Ghosh /* 349bd79ccSDebojyoti Ghosh Defines the basic matrix operations for the KAIJ matrix storage format. 449bd79ccSDebojyoti Ghosh This format is used to evaluate matrices of the form: 549bd79ccSDebojyoti Ghosh 649bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 749bd79ccSDebojyoti Ghosh 849bd79ccSDebojyoti Ghosh where 949bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 1049bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 1149bd79ccSDebojyoti Ghosh A is an AIJ (n \times n) matrix 1249bd79ccSDebojyoti Ghosh I is the identity matrix 1349bd79ccSDebojyoti Ghosh 1449bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 1549bd79ccSDebojyoti Ghosh 1649bd79ccSDebojyoti Ghosh We provide: 1749bd79ccSDebojyoti Ghosh MatMult() 1849bd79ccSDebojyoti Ghosh MatMultAdd() 1949bd79ccSDebojyoti Ghosh MatInvertBlockDiagonal() 2049bd79ccSDebojyoti Ghosh and 2149bd79ccSDebojyoti Ghosh MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*) 2249bd79ccSDebojyoti Ghosh 2349bd79ccSDebojyoti Ghosh This single directory handles both the sequential and parallel codes 2449bd79ccSDebojyoti Ghosh */ 2549bd79ccSDebojyoti Ghosh 2649bd79ccSDebojyoti Ghosh #include <../src/mat/impls/kaij/kaij.h> /*I "petscmat.h" I*/ 2749bd79ccSDebojyoti Ghosh #include <../src/mat/utils/freespace.h> 2849bd79ccSDebojyoti Ghosh #include <petsc/private/vecimpl.h> 2949bd79ccSDebojyoti Ghosh 3049bd79ccSDebojyoti Ghosh /*@C 3111a5261eSBarry Smith MatKAIJGetAIJ - Get the `MATAIJ` matrix describing the blockwise action of the `MATKAIJ` matrix 3249bd79ccSDebojyoti Ghosh 3311a5261eSBarry Smith Not Collective, but if the `MATKAIJ` matrix is parallel, the `MATAIJ` matrix is also parallel 3449bd79ccSDebojyoti Ghosh 3549bd79ccSDebojyoti Ghosh Input Parameter: 3611a5261eSBarry Smith . A - the `MATKAIJ` matrix 3749bd79ccSDebojyoti Ghosh 3849bd79ccSDebojyoti Ghosh Output Parameter: 3911a5261eSBarry Smith . B - the `MATAIJ` matrix 4049bd79ccSDebojyoti Ghosh 4149bd79ccSDebojyoti Ghosh Level: advanced 4249bd79ccSDebojyoti Ghosh 4311a5261eSBarry Smith Note: 4411a5261eSBarry Smith The reference count on the `MATAIJ` matrix is not increased so you should not destroy it. 4549bd79ccSDebojyoti Ghosh 46*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatCreateKAIJ()`, `MATKAIJ`, `MATAIJ` 4749bd79ccSDebojyoti Ghosh @*/ 48d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetAIJ(Mat A, Mat *B) 49d71ae5a4SJacob Faibussowitsch { 5049bd79ccSDebojyoti Ghosh PetscBool ismpikaij, isseqkaij; 5149bd79ccSDebojyoti Ghosh 5249bd79ccSDebojyoti Ghosh PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij)); 549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQKAIJ, &isseqkaij)); 5549bd79ccSDebojyoti Ghosh if (ismpikaij) { 5649bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 5749bd79ccSDebojyoti Ghosh 5849bd79ccSDebojyoti Ghosh *B = b->A; 5949bd79ccSDebojyoti Ghosh } else if (isseqkaij) { 6049bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 6149bd79ccSDebojyoti Ghosh 6249bd79ccSDebojyoti Ghosh *B = b->AIJ; 63b04351cbSRichard Tran Mills } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix passed in is not of type KAIJ"); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6549bd79ccSDebojyoti Ghosh } 6649bd79ccSDebojyoti Ghosh 6749bd79ccSDebojyoti Ghosh /*@C 68*2ef1f0ffSBarry Smith MatKAIJGetS - Get the `S` matrix describing the shift action of the `MATKAIJ` matrix 6949bd79ccSDebojyoti Ghosh 70*2ef1f0ffSBarry Smith Not Collective; the entire `S` is stored and returned independently on all processes. 7149bd79ccSDebojyoti Ghosh 7249bd79ccSDebojyoti Ghosh Input Parameter: 7311a5261eSBarry Smith . A - the `MATKAIJ` matrix 7449bd79ccSDebojyoti Ghosh 75a5b5c723SRichard Tran Mills Output Parameters: 76*2ef1f0ffSBarry Smith + m - the number of rows in `S` 77*2ef1f0ffSBarry Smith . n - the number of columns in `S` 78a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 7949bd79ccSDebojyoti Ghosh 8049bd79ccSDebojyoti Ghosh Level: advanced 8149bd79ccSDebojyoti Ghosh 82*2ef1f0ffSBarry Smith Note: 83*2ef1f0ffSBarry Smith All output parameters are optional (pass `NULL` if not desired) 84*2ef1f0ffSBarry Smith 85*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 8649bd79ccSDebojyoti Ghosh @*/ 87d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetS(Mat A, PetscInt *m, PetscInt *n, PetscScalar **S) 88d71ae5a4SJacob Faibussowitsch { 8949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 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 98*2ef1f0ffSBarry Smith MatKAIJGetSRead - Get a read-only pointer to the `S` matrix describing the shift action of the `MATKAIJ` matrix 99a5b5c723SRichard Tran Mills 100*2ef1f0ffSBarry 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: 106*2ef1f0ffSBarry Smith + m - the number of rows in `S` 107*2ef1f0ffSBarry 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 112*2ef1f0ffSBarry Smith Note: 113*2ef1f0ffSBarry Smith All output parameters are optional (pass `NULL` if not desired) 114*2ef1f0ffSBarry Smith 115*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 116a5b5c723SRichard Tran Mills @*/ 117d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetSRead(Mat A, PetscInt *m, PetscInt *n, const PetscScalar **S) 118d71ae5a4SJacob Faibussowitsch { 119a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 120a5b5c723SRichard Tran Mills PetscFunctionBegin; 121a5b5c723SRichard Tran Mills if (m) *m = b->p; 122a5b5c723SRichard Tran Mills if (n) *n = b->q; 123a5b5c723SRichard Tran Mills if (S) *S = b->S; 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125a5b5c723SRichard Tran Mills } 126a5b5c723SRichard Tran Mills 127a5b5c723SRichard Tran Mills /*@C 12811a5261eSBarry Smith MatKAIJRestoreS - Restore array obtained with `MatKAIJGetS()` 129a5b5c723SRichard Tran Mills 130*2ef1f0ffSBarry Smith Not Collective 131a5b5c723SRichard Tran Mills 132*2ef1f0ffSBarry Smith Input Parameters: 133*2ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix 134*2ef1f0ffSBarry Smith - S - location of pointer to array obtained with `MatKAIJGetS()` 135a5b5c723SRichard Tran Mills 136a5b5c723SRichard Tran Mills Level: advanced 13711a5261eSBarry Smith 138*2ef1f0ffSBarry Smith Note: 139*2ef1f0ffSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 140*2ef1f0ffSBarry Smith If `NULL` is passed, it will not attempt to zero the array pointer. 141*2ef1f0ffSBarry Smith 142*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`, `MatKAIJRestoreSRead()` 143a5b5c723SRichard Tran Mills @*/ 144d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreS(Mat A, PetscScalar **S) 145d71ae5a4SJacob Faibussowitsch { 146a5b5c723SRichard Tran Mills PetscFunctionBegin; 147a5b5c723SRichard Tran Mills if (S) *S = NULL; 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)A)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150a5b5c723SRichard Tran Mills } 151a5b5c723SRichard Tran Mills 152a5b5c723SRichard Tran Mills /*@C 15311a5261eSBarry Smith MatKAIJRestoreSRead - Restore array obtained with `MatKAIJGetSRead()` 154a5b5c723SRichard Tran Mills 155*2ef1f0ffSBarry Smith Not Collective 156a5b5c723SRichard Tran Mills 157*2ef1f0ffSBarry Smith Input Parameters: 158*2ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix 159*2ef1f0ffSBarry Smith - S - location of pointer to array obtained with `MatKAIJGetS()` 160a5b5c723SRichard Tran Mills 161a5b5c723SRichard Tran Mills Level: advanced 16211a5261eSBarry Smith 163*2ef1f0ffSBarry Smith Note: 164*2ef1f0ffSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 165*2ef1f0ffSBarry Smith If `NULL` is passed, it will not attempt to zero the array pointer. 166*2ef1f0ffSBarry Smith 167*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`, `MatKAIJRestoreSRead()` 168a5b5c723SRichard Tran Mills @*/ 169d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreSRead(Mat A, const PetscScalar **S) 170d71ae5a4SJacob Faibussowitsch { 171a5b5c723SRichard Tran Mills PetscFunctionBegin; 172a5b5c723SRichard Tran Mills if (S) *S = NULL; 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17449bd79ccSDebojyoti Ghosh } 17549bd79ccSDebojyoti Ghosh 17649bd79ccSDebojyoti Ghosh /*@C 177*2ef1f0ffSBarry Smith MatKAIJGetT - Get the transformation matrix `T` associated with the `MATKAIJ` matrix 17849bd79ccSDebojyoti Ghosh 179*2ef1f0ffSBarry Smith Not Collective; the entire `T` is stored and returned independently on all processes 18049bd79ccSDebojyoti Ghosh 18149bd79ccSDebojyoti Ghosh Input Parameter: 18211a5261eSBarry Smith . A - the `MATKAIJ` matrix 18349bd79ccSDebojyoti Ghosh 184d8d19677SJose E. Roman Output Parameters: 185*2ef1f0ffSBarry Smith + m - the number of rows in `T` 186*2ef1f0ffSBarry Smith . n - the number of columns in `T` 187a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 18849bd79ccSDebojyoti Ghosh 18949bd79ccSDebojyoti Ghosh Level: advanced 19049bd79ccSDebojyoti Ghosh 191*2ef1f0ffSBarry Smith Note: 192*2ef1f0ffSBarry Smith All output parameters are optional (pass `NULL` if not desired) 193*2ef1f0ffSBarry Smith 194*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 19549bd79ccSDebojyoti Ghosh @*/ 196d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetT(Mat A, PetscInt *m, PetscInt *n, PetscScalar **T) 197d71ae5a4SJacob Faibussowitsch { 19849bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 19949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 200a5b5c723SRichard Tran Mills if (m) *m = b->p; 201a5b5c723SRichard Tran Mills if (n) *n = b->q; 202a5b5c723SRichard Tran Mills if (T) *T = b->T; 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204a5b5c723SRichard Tran Mills } 205a5b5c723SRichard Tran Mills 206a5b5c723SRichard Tran Mills /*@C 207*2ef1f0ffSBarry Smith MatKAIJGetTRead - Get a read-only pointer to the transformation matrix `T` associated with the `MATKAIJ` matrix 208a5b5c723SRichard Tran Mills 209*2ef1f0ffSBarry Smith Not Collective; the entire `T` is stored and returned independently on all processes 210a5b5c723SRichard Tran Mills 211a5b5c723SRichard Tran Mills Input Parameter: 21211a5261eSBarry Smith . A - the `MATKAIJ` matrix 213a5b5c723SRichard Tran Mills 214d8d19677SJose E. Roman Output Parameters: 215*2ef1f0ffSBarry Smith + m - the number of rows in `T` 216*2ef1f0ffSBarry Smith . n - the number of columns in `T` 217a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 218a5b5c723SRichard Tran Mills 219a5b5c723SRichard Tran Mills Level: advanced 220a5b5c723SRichard Tran Mills 221*2ef1f0ffSBarry Smith Note: 222*2ef1f0ffSBarry Smith All output parameters are optional (pass `NULL` if not desired) 223*2ef1f0ffSBarry Smith 224*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 225a5b5c723SRichard Tran Mills @*/ 226d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetTRead(Mat A, PetscInt *m, PetscInt *n, const PetscScalar **T) 227d71ae5a4SJacob Faibussowitsch { 228a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 229a5b5c723SRichard Tran Mills PetscFunctionBegin; 230a5b5c723SRichard Tran Mills if (m) *m = b->p; 231a5b5c723SRichard Tran Mills if (n) *n = b->q; 232a5b5c723SRichard Tran Mills if (T) *T = b->T; 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 234a5b5c723SRichard Tran Mills } 235a5b5c723SRichard Tran Mills 236a5b5c723SRichard Tran Mills /*@C 23711a5261eSBarry Smith MatKAIJRestoreT - Restore array obtained with `MatKAIJGetT()` 238a5b5c723SRichard Tran Mills 239*2ef1f0ffSBarry Smith Not Collective 240a5b5c723SRichard Tran Mills 241*2ef1f0ffSBarry Smith Input Parameters: 242*2ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix 243*2ef1f0ffSBarry Smith - T - location of pointer to array obtained with `MatKAIJGetS()` 244a5b5c723SRichard Tran Mills 245a5b5c723SRichard Tran Mills Level: advanced 24611a5261eSBarry Smith 247*2ef1f0ffSBarry Smith Note: 248*2ef1f0ffSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 249*2ef1f0ffSBarry Smith If `NULL` is passed, it will not attempt to zero the array pointer. 250*2ef1f0ffSBarry Smith 251*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`, `MatKAIJRestoreTRead()` 252a5b5c723SRichard Tran Mills @*/ 253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreT(Mat A, PetscScalar **T) 254d71ae5a4SJacob Faibussowitsch { 255a5b5c723SRichard Tran Mills PetscFunctionBegin; 256a5b5c723SRichard Tran Mills if (T) *T = NULL; 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)A)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259a5b5c723SRichard Tran Mills } 260a5b5c723SRichard Tran Mills 261a5b5c723SRichard Tran Mills /*@C 26211a5261eSBarry Smith MatKAIJRestoreTRead - Restore array obtained with `MatKAIJGetTRead()` 263a5b5c723SRichard Tran Mills 264*2ef1f0ffSBarry Smith Not Collective 265a5b5c723SRichard Tran Mills 266*2ef1f0ffSBarry Smith Input Parameters: 267*2ef1f0ffSBarry Smith + A - the `MATKAIJ` matrix 268*2ef1f0ffSBarry Smith - T - location of pointer to array obtained with `MatKAIJGetS()` 269a5b5c723SRichard Tran Mills 270a5b5c723SRichard Tran Mills Level: advanced 27111a5261eSBarry Smith 272*2ef1f0ffSBarry Smith Note: 273*2ef1f0ffSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 274*2ef1f0ffSBarry Smith If `NULL` is passed, it will not attempt to zero the array pointer. 275*2ef1f0ffSBarry Smith 276*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`, `MatKAIJRestoreTRead()` 277a5b5c723SRichard Tran Mills @*/ 278d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreTRead(Mat A, const PetscScalar **T) 279d71ae5a4SJacob Faibussowitsch { 280a5b5c723SRichard Tran Mills PetscFunctionBegin; 281a5b5c723SRichard Tran Mills if (T) *T = NULL; 2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28349bd79ccSDebojyoti Ghosh } 28449bd79ccSDebojyoti Ghosh 2850567c835SRichard Tran Mills /*@ 28611a5261eSBarry Smith MatKAIJSetAIJ - Set the `MATAIJ` matrix describing the blockwise action of the `MATKAIJ` matrix 2870567c835SRichard Tran Mills 28811a5261eSBarry Smith Logically Collective; if the `MATAIJ` matrix is parallel, the `MATKAIJ` matrix is also parallel 2890567c835SRichard Tran Mills 2900567c835SRichard Tran Mills Input Parameters: 29111a5261eSBarry Smith + A - the `MATKAIJ` matrix 29211a5261eSBarry Smith - B - the `MATAIJ` matrix 2930567c835SRichard Tran Mills 294*2ef1f0ffSBarry Smith Level: advanced 295*2ef1f0ffSBarry Smith 29615b9d025SRichard Tran Mills Notes: 29711a5261eSBarry 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. 29811a5261eSBarry Smith 29911a5261eSBarry Smith Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix. 30015b9d025SRichard Tran Mills 301*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()` 3020567c835SRichard Tran Mills @*/ 303d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetAIJ(Mat A, Mat B) 304d71ae5a4SJacob Faibussowitsch { 3050567c835SRichard Tran Mills PetscMPIInt size; 30606d61a37SPierre Jolivet PetscBool flg; 3070567c835SRichard Tran Mills 3080567c835SRichard Tran Mills PetscFunctionBegin; 3099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 3100567c835SRichard Tran Mills if (size == 1) { 3119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg)); 31228b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MatKAIJSetAIJ() with MATSEQKAIJ does not support %s as the AIJ mat", ((PetscObject)B)->type_name); 3130567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 3140567c835SRichard Tran Mills a->AIJ = B; 3150567c835SRichard Tran Mills } else { 3160567c835SRichard Tran Mills Mat_MPIKAIJ *a = (Mat_MPIKAIJ *)A->data; 3170567c835SRichard Tran Mills a->A = B; 3180567c835SRichard Tran Mills } 3199566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3210567c835SRichard Tran Mills } 3220567c835SRichard Tran Mills 3230567c835SRichard Tran Mills /*@C 324*2ef1f0ffSBarry Smith MatKAIJSetS - Set the `S` matrix describing the shift action of the `MATKAIJ` matrix 3250567c835SRichard Tran Mills 326*2ef1f0ffSBarry Smith Logically Collective; the entire `S` is stored independently on all processes. 3270567c835SRichard Tran Mills 3280567c835SRichard Tran Mills Input Parameters: 32911a5261eSBarry Smith + A - the `MATKAIJ` matrix 330*2ef1f0ffSBarry Smith . p - the number of rows in `S` 331*2ef1f0ffSBarry Smith . q - the number of columns in `S` 3320567c835SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 3330567c835SRichard Tran Mills 33411a5261eSBarry Smith Level: advanced 3350567c835SRichard Tran Mills 336*2ef1f0ffSBarry Smith Notes: 337*2ef1f0ffSBarry Smith The dimensions `p` and `q` must match those of the transformation matrix `T` associated with the `MATKAIJ` matrix. 338*2ef1f0ffSBarry Smith 339*2ef1f0ffSBarry Smith The `S` matrix is copied, so the user can destroy this array. 340*2ef1f0ffSBarry Smith 341*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJSetT()`, `MatKAIJSetAIJ()` 3420567c835SRichard Tran Mills @*/ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetS(Mat A, PetscInt p, PetscInt q, const PetscScalar S[]) 344d71ae5a4SJacob Faibussowitsch { 3450567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 3460567c835SRichard Tran Mills 3470567c835SRichard Tran Mills PetscFunctionBegin; 3489566063dSJacob Faibussowitsch PetscCall(PetscFree(a->S)); 3490567c835SRichard Tran Mills if (S) { 3509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p * q, &a->S)); 3519566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(a->S, S, p * q * sizeof(PetscScalar))); 3520567c835SRichard Tran Mills } else a->S = NULL; 3530567c835SRichard Tran Mills 3540567c835SRichard Tran Mills a->p = p; 3550567c835SRichard Tran Mills a->q = q; 3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3570567c835SRichard Tran Mills } 3580567c835SRichard Tran Mills 3590567c835SRichard Tran Mills /*@C 360*2ef1f0ffSBarry Smith MatKAIJGetScaledIdentity - Check if both `S` and `T` are scaled identities. 361910cf402Sprj- 362910cf402Sprj- Logically Collective. 363910cf402Sprj- 364910cf402Sprj- Input Parameter: 36511a5261eSBarry Smith . A - the `MATKAIJ` matrix 366910cf402Sprj- 367910cf402Sprj- Output Parameter: 368910cf402Sprj- . identity - the Boolean value 369910cf402Sprj- 370910cf402Sprj- Level: Advanced 371910cf402Sprj- 372*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetT()` 373910cf402Sprj- @*/ 374d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetScaledIdentity(Mat A, PetscBool *identity) 375d71ae5a4SJacob Faibussowitsch { 376910cf402Sprj- Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 377910cf402Sprj- PetscInt i, j; 378910cf402Sprj- 379910cf402Sprj- PetscFunctionBegin; 380910cf402Sprj- if (a->p != a->q) { 381910cf402Sprj- *identity = PETSC_FALSE; 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383910cf402Sprj- } else *identity = PETSC_TRUE; 384910cf402Sprj- if (!a->isTI || a->S) { 385910cf402Sprj- for (i = 0; i < a->p && *identity; i++) { 386910cf402Sprj- for (j = 0; j < a->p && *identity; j++) { 387910cf402Sprj- if (i != j) { 388910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 389910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 390910cf402Sprj- } else { 391910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i * (a->p + 1)] - a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 392910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i * (a->p + 1)] - a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 393910cf402Sprj- } 394910cf402Sprj- } 395910cf402Sprj- } 396910cf402Sprj- } 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 398910cf402Sprj- } 399910cf402Sprj- 400910cf402Sprj- /*@C 401*2ef1f0ffSBarry Smith MatKAIJSetT - Set the transformation matrix `T` associated with the `MATKAIJ` matrix 4020567c835SRichard Tran Mills 403*2ef1f0ffSBarry Smith Logically Collective; the entire `T` is stored independently on all processes. 4040567c835SRichard Tran Mills 4050567c835SRichard Tran Mills Input Parameters: 40611a5261eSBarry Smith + A - the `MATKAIJ` matrix 407*2ef1f0ffSBarry Smith . p - the number of rows in `S` 408*2ef1f0ffSBarry Smith . q - the number of columns in `S` 409*2ef1f0ffSBarry Smith - T - the `T` matrix, in form of a scalar array in column-major format 4100567c835SRichard Tran Mills 4110567c835SRichard Tran Mills Level: Advanced 4120567c835SRichard Tran Mills 413*2ef1f0ffSBarry Smith Notes: 414*2ef1f0ffSBarry Smith The dimensions `p` and `q` must match those of the shift matrix `S` associated with the `MATKAIJ` matrix. 415*2ef1f0ffSBarry Smith 416*2ef1f0ffSBarry Smith The `T` matrix is copied, so the user can destroy this array. 417*2ef1f0ffSBarry Smith 418*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJSetS()`, `MatKAIJSetAIJ()` 4190567c835SRichard Tran Mills @*/ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetT(Mat A, PetscInt p, PetscInt q, const PetscScalar T[]) 421d71ae5a4SJacob Faibussowitsch { 4220567c835SRichard Tran Mills PetscInt i, j; 4230567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 4240567c835SRichard Tran Mills PetscBool isTI = PETSC_FALSE; 4250567c835SRichard Tran Mills 4260567c835SRichard Tran Mills PetscFunctionBegin; 4270567c835SRichard Tran Mills /* check if T is an identity matrix */ 4280567c835SRichard Tran Mills if (T && (p == q)) { 4290567c835SRichard Tran Mills isTI = PETSC_TRUE; 4300567c835SRichard Tran Mills for (i = 0; i < p; i++) { 4310567c835SRichard Tran Mills for (j = 0; j < q; j++) { 4320567c835SRichard Tran Mills if (i == j) { 4330567c835SRichard Tran Mills /* diagonal term must be 1 */ 4340567c835SRichard Tran Mills if (T[i + j * p] != 1.0) isTI = PETSC_FALSE; 4350567c835SRichard Tran Mills } else { 4360567c835SRichard Tran Mills /* off-diagonal term must be 0 */ 4370567c835SRichard Tran Mills if (T[i + j * p] != 0.0) isTI = PETSC_FALSE; 4380567c835SRichard Tran Mills } 4390567c835SRichard Tran Mills } 4400567c835SRichard Tran Mills } 4410567c835SRichard Tran Mills } 4420567c835SRichard Tran Mills a->isTI = isTI; 4430567c835SRichard Tran Mills 4449566063dSJacob Faibussowitsch PetscCall(PetscFree(a->T)); 4450567c835SRichard Tran Mills if (T && (!isTI)) { 4469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p * q, &a->T)); 4479566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(a->T, T, p * q * sizeof(PetscScalar))); 44850d19d74SRichard Tran Mills } else a->T = NULL; 4490567c835SRichard Tran Mills 4500567c835SRichard Tran Mills a->p = p; 4510567c835SRichard Tran Mills a->q = q; 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4530567c835SRichard Tran Mills } 4540567c835SRichard Tran Mills 455ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_SeqKAIJ(Mat A) 456d71ae5a4SJacob Faibussowitsch { 45749bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 45849bd79ccSDebojyoti Ghosh 45949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 4609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->AIJ)); 4619566063dSJacob Faibussowitsch PetscCall(PetscFree(b->S)); 4629566063dSJacob Faibussowitsch PetscCall(PetscFree(b->T)); 4639566063dSJacob Faibussowitsch PetscCall(PetscFree(b->ibdiag)); 4649566063dSJacob Faibussowitsch PetscCall(PetscFree5(b->sor.w, b->sor.y, b->sor.work, b->sor.t, b->sor.arr)); 4659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", NULL)); 4669566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46849bd79ccSDebojyoti Ghosh } 46949bd79ccSDebojyoti Ghosh 470ff6a9541SJacob Faibussowitsch static PetscErrorCode MatKAIJ_build_AIJ_OAIJ(Mat A) 471d71ae5a4SJacob Faibussowitsch { 472e0e5a793SRichard Tran Mills Mat_MPIKAIJ *a; 473e0e5a793SRichard Tran Mills Mat_MPIAIJ *mpiaij; 474e0e5a793SRichard Tran Mills PetscScalar *T; 475e0e5a793SRichard Tran Mills PetscInt i, j; 476e0e5a793SRichard Tran Mills PetscObjectState state; 477e0e5a793SRichard Tran Mills 478e0e5a793SRichard Tran Mills PetscFunctionBegin; 479e0e5a793SRichard Tran Mills a = (Mat_MPIKAIJ *)A->data; 480e0e5a793SRichard Tran Mills mpiaij = (Mat_MPIAIJ *)a->A->data; 481e0e5a793SRichard Tran Mills 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)a->A, &state)); 483e0e5a793SRichard Tran Mills if (state == a->state) { 484e0e5a793SRichard Tran Mills /* The existing AIJ and KAIJ members are up-to-date, so simply exit. */ 4853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 486e0e5a793SRichard Tran Mills } else { 4879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&a->AIJ)); 4889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&a->OAIJ)); 489e0e5a793SRichard Tran Mills if (a->isTI) { 490e0e5a793SRichard Tran Mills /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL. 491e0e5a793SRichard Tran Mills * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will 492e0e5a793SRichard Tran Mills * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix 493e0e5a793SRichard Tran Mills * to pass in. */ 4949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->p * a->q, &T)); 495e0e5a793SRichard Tran Mills for (i = 0; i < a->p; i++) { 496e0e5a793SRichard Tran Mills for (j = 0; j < a->q; j++) { 497e0e5a793SRichard Tran Mills if (i == j) T[i + j * a->p] = 1.0; 498e0e5a793SRichard Tran Mills else T[i + j * a->p] = 0.0; 499e0e5a793SRichard Tran Mills } 500e0e5a793SRichard Tran Mills } 501e0e5a793SRichard Tran Mills } else T = a->T; 5029566063dSJacob Faibussowitsch PetscCall(MatCreateKAIJ(mpiaij->A, a->p, a->q, a->S, T, &a->AIJ)); 5039566063dSJacob Faibussowitsch PetscCall(MatCreateKAIJ(mpiaij->B, a->p, a->q, NULL, T, &a->OAIJ)); 5041baa6e33SBarry Smith if (a->isTI) PetscCall(PetscFree(T)); 505e0e5a793SRichard Tran Mills a->state = state; 506e0e5a793SRichard Tran Mills } 507e0e5a793SRichard Tran Mills 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509e0e5a793SRichard Tran Mills } 510e0e5a793SRichard Tran Mills 511ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSetUp_KAIJ(Mat A) 512d71ae5a4SJacob Faibussowitsch { 5130567c835SRichard Tran Mills PetscInt n; 5140567c835SRichard Tran Mills PetscMPIInt size; 5150567c835SRichard Tran Mills Mat_SeqKAIJ *seqkaij = (Mat_SeqKAIJ *)A->data; 5160567c835SRichard Tran Mills 51749bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 5190567c835SRichard Tran Mills if (size == 1) { 5209566063dSJacob 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)); 5219566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p)); 5229566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q)); 5239566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5249566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5250567c835SRichard Tran Mills } else { 5260567c835SRichard Tran Mills Mat_MPIKAIJ *a; 5270567c835SRichard Tran Mills Mat_MPIAIJ *mpiaij; 5280567c835SRichard Tran Mills IS from, to; 5290567c835SRichard Tran Mills Vec gvec; 5300567c835SRichard Tran Mills 5310567c835SRichard Tran Mills a = (Mat_MPIKAIJ *)A->data; 532d3f912faSRichard Tran Mills mpiaij = (Mat_MPIAIJ *)a->A->data; 5339566063dSJacob 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)); 5349566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p)); 5359566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q)); 5369566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5379566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5380567c835SRichard Tran Mills 5399566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); 5400567c835SRichard Tran Mills 5419566063dSJacob Faibussowitsch PetscCall(VecGetSize(mpiaij->lvec, &n)); 5429566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &a->w)); 5439566063dSJacob Faibussowitsch PetscCall(VecSetSizes(a->w, n * a->q, n * a->q)); 5449566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(a->w, a->q)); 5459566063dSJacob Faibussowitsch PetscCall(VecSetType(a->w, VECSEQ)); 5460567c835SRichard Tran Mills 5470567c835SRichard Tran Mills /* create two temporary Index sets for build scatter gather */ 5489566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)a->A), a->q, n, mpiaij->garray, PETSC_COPY_VALUES, &from)); 5499566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, n * a->q, 0, 1, &to)); 5500567c835SRichard Tran Mills 5510567c835SRichard Tran Mills /* create temporary global vector to generate scatter context */ 5529566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A), a->q, a->q * a->A->cmap->n, a->q * a->A->cmap->N, NULL, &gvec)); 5530567c835SRichard Tran Mills 5540567c835SRichard Tran Mills /* generate the scatter context */ 5559566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(gvec, from, a->w, to, &a->ctx)); 5560567c835SRichard Tran Mills 5579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 5589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 5599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gvec)); 5600567c835SRichard Tran Mills } 5610567c835SRichard Tran Mills 5620567c835SRichard Tran Mills A->assembled = PETSC_TRUE; 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56449bd79ccSDebojyoti Ghosh } 56549bd79ccSDebojyoti Ghosh 566ff6a9541SJacob Faibussowitsch static PetscErrorCode MatView_KAIJ(Mat A, PetscViewer viewer) 567d71ae5a4SJacob Faibussowitsch { 568e6985dafSRichard Tran Mills PetscViewerFormat format; 569e6985dafSRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 57049bd79ccSDebojyoti Ghosh Mat B; 571e6985dafSRichard Tran Mills PetscInt i; 572e6985dafSRichard Tran Mills PetscBool ismpikaij; 57349bd79ccSDebojyoti Ghosh 57449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij)); 5769566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 577e6985dafSRichard Tran Mills if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) { 5789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "S and T have %" PetscInt_FMT " rows and %" PetscInt_FMT " columns\n", a->p, a->q)); 579e6985dafSRichard Tran Mills 580e6985dafSRichard Tran Mills /* Print appropriate details for S. */ 581e6985dafSRichard Tran Mills if (!a->S) { 5829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "S is NULL\n")); 583e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 5849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of S are ")); 585e6985dafSRichard Tran Mills for (i = 0; i < (a->p * a->q); i++) { 586e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 5879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->S[i]), (double)PetscImaginaryPart(a->S[i]))); 588e6985dafSRichard Tran Mills #else 5899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->S[i]))); 590e6985dafSRichard Tran Mills #endif 591e6985dafSRichard Tran Mills } 5929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 59349bd79ccSDebojyoti Ghosh } 59449bd79ccSDebojyoti Ghosh 595e6985dafSRichard Tran Mills /* Print appropriate details for T. */ 596e6985dafSRichard Tran Mills if (a->isTI) { 5979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "T is the identity matrix\n")); 598e6985dafSRichard Tran Mills } else if (!a->T) { 5999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "T is NULL\n")); 600e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 6019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of T are ")); 602e6985dafSRichard Tran Mills for (i = 0; i < (a->p * a->q); i++) { 603e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 6049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->T[i]), (double)PetscImaginaryPart(a->T[i]))); 605e6985dafSRichard Tran Mills #else 6069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->T[i]))); 607e6985dafSRichard Tran Mills #endif 608e6985dafSRichard Tran Mills } 6099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 610e6985dafSRichard Tran Mills } 61149bd79ccSDebojyoti Ghosh 612e6985dafSRichard Tran Mills /* Now print details for the AIJ matrix, using the AIJ viewer. */ 6139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Now viewing the associated AIJ matrix:\n")); 614e6985dafSRichard Tran Mills if (ismpikaij) { 615e6985dafSRichard Tran Mills Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 6169566063dSJacob Faibussowitsch PetscCall(MatView(b->A, viewer)); 617e6985dafSRichard Tran Mills } else { 6189566063dSJacob Faibussowitsch PetscCall(MatView(a->AIJ, viewer)); 619e6985dafSRichard Tran Mills } 620e6985dafSRichard Tran Mills 621e6985dafSRichard Tran Mills } else { 622e6985dafSRichard Tran Mills /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */ 6239566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B)); 6249566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 6259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 626e6985dafSRichard Tran Mills } 6273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62849bd79ccSDebojyoti Ghosh } 62949bd79ccSDebojyoti Ghosh 630ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_MPIKAIJ(Mat A) 631d71ae5a4SJacob Faibussowitsch { 63249bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 63349bd79ccSDebojyoti Ghosh 63449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 6359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->AIJ)); 6369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->OAIJ)); 6379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->A)); 6389566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&b->ctx)); 6399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->w)); 6409566063dSJacob Faibussowitsch PetscCall(PetscFree(b->S)); 6419566063dSJacob Faibussowitsch PetscCall(PetscFree(b->T)); 6429566063dSJacob Faibussowitsch PetscCall(PetscFree(b->ibdiag)); 6439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", NULL)); 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", NULL)); 6459566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 6463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64749bd79ccSDebojyoti Ghosh } 64849bd79ccSDebojyoti Ghosh 64949bd79ccSDebojyoti Ghosh /* zz = yy + Axx */ 650ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_SeqKAIJ(Mat A, Vec xx, Vec yy, Vec zz) 651d71ae5a4SJacob Faibussowitsch { 65249bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 65349bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; 65449bd79ccSDebojyoti Ghosh const PetscScalar *s = b->S, *t = b->T; 65549bd79ccSDebojyoti Ghosh const PetscScalar *x, *v, *bx; 65649bd79ccSDebojyoti Ghosh PetscScalar *y, *sums; 65749bd79ccSDebojyoti Ghosh const PetscInt m = b->AIJ->rmap->n, *idx, *ii; 65849bd79ccSDebojyoti Ghosh PetscInt n, i, jrow, j, l, p = b->p, q = b->q, k; 65949bd79ccSDebojyoti Ghosh 66049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 66149bd79ccSDebojyoti Ghosh if (!yy) { 6629566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 66349bd79ccSDebojyoti Ghosh } else { 6649566063dSJacob Faibussowitsch PetscCall(VecCopy(yy, zz)); 66549bd79ccSDebojyoti Ghosh } 6663ba16761SJacob Faibussowitsch if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(PETSC_SUCCESS); 66749bd79ccSDebojyoti Ghosh 6689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x)); 6699566063dSJacob Faibussowitsch PetscCall(VecGetArray(zz, &y)); 67049bd79ccSDebojyoti Ghosh idx = a->j; 67149bd79ccSDebojyoti Ghosh v = a->a; 67249bd79ccSDebojyoti Ghosh ii = a->i; 67349bd79ccSDebojyoti Ghosh 67449bd79ccSDebojyoti Ghosh if (b->isTI) { 67549bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 67649bd79ccSDebojyoti Ghosh jrow = ii[i]; 67749bd79ccSDebojyoti Ghosh n = ii[i + 1] - jrow; 67849bd79ccSDebojyoti Ghosh sums = y + p * i; 67949bd79ccSDebojyoti Ghosh for (j = 0; j < n; j++) { 680ad540459SPierre Jolivet for (k = 0; k < p; k++) sums[k] += v[jrow + j] * x[q * idx[jrow + j] + k]; 68149bd79ccSDebojyoti Ghosh } 68249bd79ccSDebojyoti Ghosh } 6839566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(3.0 * (a->nz) * p)); 68449bd79ccSDebojyoti Ghosh } else if (t) { 68549bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 68649bd79ccSDebojyoti Ghosh jrow = ii[i]; 68749bd79ccSDebojyoti Ghosh n = ii[i + 1] - jrow; 68849bd79ccSDebojyoti Ghosh sums = y + p * i; 68949bd79ccSDebojyoti Ghosh for (j = 0; j < n; j++) { 69049bd79ccSDebojyoti Ghosh for (k = 0; k < p; k++) { 691ad540459SPierre Jolivet for (l = 0; l < q; l++) sums[k] += v[jrow + j] * t[k + l * p] * x[q * idx[jrow + j] + l]; 69249bd79ccSDebojyoti Ghosh } 69349bd79ccSDebojyoti Ghosh } 69449bd79ccSDebojyoti Ghosh } 695234d9204SRichard Tran Mills /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do), 696234d9204SRichard Tran Mills * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T), 697234d9204SRichard Tran Mills * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter 698234d9204SRichard Tran Mills * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */ 6999566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((2.0 * p * q - p) * m + 2.0 * p * a->nz)); 70049bd79ccSDebojyoti Ghosh } 70149bd79ccSDebojyoti Ghosh if (s) { 70249bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 70349bd79ccSDebojyoti Ghosh sums = y + p * i; 70449bd79ccSDebojyoti Ghosh bx = x + q * i; 70549bd79ccSDebojyoti Ghosh if (i < b->AIJ->cmap->n) { 70649bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 707ad540459SPierre Jolivet for (k = 0; k < p; k++) sums[k] += s[k + j * p] * bx[j]; 70849bd79ccSDebojyoti Ghosh } 70949bd79ccSDebojyoti Ghosh } 71049bd79ccSDebojyoti Ghosh } 7119566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * m * p * q)); 71249bd79ccSDebojyoti Ghosh } 71349bd79ccSDebojyoti Ghosh 7149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x)); 7159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(zz, &y)); 7163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71749bd79ccSDebojyoti Ghosh } 71849bd79ccSDebojyoti Ghosh 719ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_SeqKAIJ(Mat A, Vec xx, Vec yy) 720d71ae5a4SJacob Faibussowitsch { 72149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 722f3fa974cSJacob Faibussowitsch PetscCall(MatMultAdd_SeqKAIJ(A, xx, NULL, yy)); 7233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72449bd79ccSDebojyoti Ghosh } 72549bd79ccSDebojyoti Ghosh 72649bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h> 72749bd79ccSDebojyoti Ghosh 728ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A, const PetscScalar **values) 729d71ae5a4SJacob Faibussowitsch { 73049bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 73149bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; 73249bd79ccSDebojyoti Ghosh const PetscScalar *S = b->S; 73349bd79ccSDebojyoti Ghosh const PetscScalar *T = b->T; 73449bd79ccSDebojyoti Ghosh const PetscScalar *v = a->a; 73549bd79ccSDebojyoti Ghosh const PetscInt p = b->p, q = b->q, m = b->AIJ->rmap->n, *idx = a->j, *ii = a->i; 73649bd79ccSDebojyoti Ghosh PetscInt i, j, *v_pivots, dof, dof2; 73749bd79ccSDebojyoti Ghosh PetscScalar *diag, aval, *v_work; 73849bd79ccSDebojyoti Ghosh 73949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 74008401ef6SPierre Jolivet PetscCheck(p == q, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Block size must be square to calculate inverse."); 741aed4548fSBarry Smith PetscCheck(S || T || b->isTI, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Cannot invert a zero matrix."); 74249bd79ccSDebojyoti Ghosh 74349bd79ccSDebojyoti Ghosh dof = p; 74449bd79ccSDebojyoti Ghosh dof2 = dof * dof; 74549bd79ccSDebojyoti Ghosh 74649bd79ccSDebojyoti Ghosh if (b->ibdiagvalid) { 74749bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 7483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74949bd79ccSDebojyoti Ghosh } 750aa624791SPierre Jolivet if (!b->ibdiag) PetscCall(PetscMalloc1(dof2 * m, &b->ibdiag)); 75149bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 75249bd79ccSDebojyoti Ghosh diag = b->ibdiag; 75349bd79ccSDebojyoti Ghosh 7549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dof, &v_work, dof, &v_pivots)); 75549bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 75649bd79ccSDebojyoti Ghosh if (S) { 7579566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(diag, S, dof2 * sizeof(PetscScalar))); 75849bd79ccSDebojyoti Ghosh } else { 7599566063dSJacob Faibussowitsch PetscCall(PetscMemzero(diag, dof2 * sizeof(PetscScalar))); 76049bd79ccSDebojyoti Ghosh } 76149bd79ccSDebojyoti Ghosh if (b->isTI) { 76249bd79ccSDebojyoti Ghosh aval = 0; 7639371c9d4SSatish Balay for (j = ii[i]; j < ii[i + 1]; j++) 7649371c9d4SSatish Balay if (idx[j] == i) aval = v[j]; 76549bd79ccSDebojyoti Ghosh for (j = 0; j < dof; j++) diag[j + dof * j] += aval; 76649bd79ccSDebojyoti Ghosh } else if (T) { 76749bd79ccSDebojyoti Ghosh aval = 0; 7689371c9d4SSatish Balay for (j = ii[i]; j < ii[i + 1]; j++) 7699371c9d4SSatish Balay if (idx[j] == i) aval = v[j]; 77049bd79ccSDebojyoti Ghosh for (j = 0; j < dof2; j++) diag[j] += aval * T[j]; 77149bd79ccSDebojyoti Ghosh } 7729566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A(dof, diag, v_pivots, v_work, PETSC_FALSE, NULL)); 77349bd79ccSDebojyoti Ghosh diag += dof2; 77449bd79ccSDebojyoti Ghosh } 7759566063dSJacob Faibussowitsch PetscCall(PetscFree2(v_work, v_pivots)); 77649bd79ccSDebojyoti Ghosh 77749bd79ccSDebojyoti Ghosh b->ibdiagvalid = PETSC_TRUE; 7783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77949bd79ccSDebojyoti Ghosh } 78049bd79ccSDebojyoti Ghosh 781d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A, Mat *B) 782d71ae5a4SJacob Faibussowitsch { 78349bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ *)A->data; 78449bd79ccSDebojyoti Ghosh 78549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 78649bd79ccSDebojyoti Ghosh *B = kaij->AIJ; 7873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78849bd79ccSDebojyoti Ghosh } 78949bd79ccSDebojyoti Ghosh 790d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConvert_KAIJ_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 791d71ae5a4SJacob Faibussowitsch { 792fac53fa1SPierre Jolivet Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 793fac53fa1SPierre Jolivet Mat AIJ, OAIJ, B; 794fac53fa1SPierre Jolivet PetscInt *d_nnz, *o_nnz = NULL, nz, i, j, m, d; 795fac53fa1SPierre Jolivet const PetscInt p = a->p, q = a->q; 796fac53fa1SPierre Jolivet PetscBool ismpikaij, missing; 797fac53fa1SPierre Jolivet 798fac53fa1SPierre Jolivet PetscFunctionBegin; 799fac53fa1SPierre Jolivet if (reuse != MAT_REUSE_MATRIX) { 8009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij)); 801fac53fa1SPierre Jolivet if (ismpikaij) { 802fac53fa1SPierre Jolivet Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 803fac53fa1SPierre Jolivet AIJ = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ; 804fac53fa1SPierre Jolivet OAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ; 805fac53fa1SPierre Jolivet } else { 806fac53fa1SPierre Jolivet AIJ = a->AIJ; 807fac53fa1SPierre Jolivet OAIJ = NULL; 808fac53fa1SPierre Jolivet } 8099566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 8109566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 8119566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 8129566063dSJacob Faibussowitsch PetscCall(MatGetSize(AIJ, &m, NULL)); 8139566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(AIJ, &missing, &d)); /* assumption that all successive rows will have a missing diagonal */ 814fac53fa1SPierre Jolivet if (!missing || !a->S) d = m; 8159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &d_nnz)); 816fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8179566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(AIJ, i, &nz, NULL, NULL)); 818fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) d_nnz[i * p + j] = nz * q + (i >= d) * q; 8199566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(AIJ, i, &nz, NULL, NULL)); 820fac53fa1SPierre Jolivet } 821fac53fa1SPierre Jolivet if (OAIJ) { 8229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &o_nnz)); 823fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8249566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL)); 825fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) o_nnz[i * p + j] = nz * q; 8269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL)); 827fac53fa1SPierre Jolivet } 8289566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B, 0, d_nnz, 0, o_nnz)); 829fac53fa1SPierre Jolivet } else { 8309566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(B, 0, d_nnz)); 831fac53fa1SPierre Jolivet } 8329566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz)); 8339566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nnz)); 834fac53fa1SPierre Jolivet } else B = *newmat; 8359566063dSJacob Faibussowitsch PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B)); 836fac53fa1SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 8379566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 838fac53fa1SPierre Jolivet } else *newmat = B; 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 840fac53fa1SPierre Jolivet } 841fac53fa1SPierre Jolivet 842ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSOR_SeqKAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 843d71ae5a4SJacob Faibussowitsch { 84449bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *kaij = (Mat_SeqKAIJ *)A->data; 84549bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)kaij->AIJ->data; 84649bd79ccSDebojyoti Ghosh const PetscScalar *aa = a->a, *T = kaij->T, *v; 84749bd79ccSDebojyoti Ghosh const PetscInt m = kaij->AIJ->rmap->n, *ai = a->i, *aj = a->j, p = kaij->p, q = kaij->q, *diag, *vi; 84849bd79ccSDebojyoti Ghosh const PetscScalar *b, *xb, *idiag; 84949bd79ccSDebojyoti Ghosh PetscScalar *x, *work, *workt, *w, *y, *arr, *t, *arrt; 85049bd79ccSDebojyoti Ghosh PetscInt i, j, k, i2, bs, bs2, nz; 85149bd79ccSDebojyoti Ghosh 85249bd79ccSDebojyoti Ghosh PetscFunctionBegin; 85349bd79ccSDebojyoti Ghosh its = its * lits; 854aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat"); 85508401ef6SPierre 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); 85628b400f6SJacob Faibussowitsch PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for diagonal shift"); 85708401ef6SPierre 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"); 85808401ef6SPierre Jolivet PetscCheck(p == q, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatSOR for KAIJ: No support for non-square dense blocks"); 859f7d195e4SLawrence Mitchell bs = p; 860f7d195e4SLawrence Mitchell bs2 = bs * bs; 86149bd79ccSDebojyoti Ghosh 8623ba16761SJacob Faibussowitsch if (!m) PetscFunctionReturn(PETSC_SUCCESS); 86349bd79ccSDebojyoti Ghosh 8649566063dSJacob Faibussowitsch if (!kaij->ibdiagvalid) PetscCall(MatInvertBlockDiagonal_SeqKAIJ(A, NULL)); 86549bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 86649bd79ccSDebojyoti Ghosh diag = a->diag; 86749bd79ccSDebojyoti Ghosh 86849bd79ccSDebojyoti Ghosh if (!kaij->sor.setup) { 8699566063dSJacob 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)); 87049bd79ccSDebojyoti Ghosh kaij->sor.setup = PETSC_TRUE; 87149bd79ccSDebojyoti Ghosh } 87249bd79ccSDebojyoti Ghosh y = kaij->sor.y; 87349bd79ccSDebojyoti Ghosh w = kaij->sor.w; 87449bd79ccSDebojyoti Ghosh work = kaij->sor.work; 87549bd79ccSDebojyoti Ghosh t = kaij->sor.t; 87649bd79ccSDebojyoti Ghosh arr = kaij->sor.arr; 87749bd79ccSDebojyoti Ghosh 878d0609cedSBarry Smith PetscCall(VecGetArray(xx, &x)); 8799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 88049bd79ccSDebojyoti Ghosh 88149bd79ccSDebojyoti Ghosh if (flag & SOR_ZERO_INITIAL_GUESS) { 88249bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 88349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, b, idiag, x); /* x[0:bs] <- D^{-1} b[0:bs] */ 8849566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t, b, bs * sizeof(PetscScalar))); 88549bd79ccSDebojyoti Ghosh i2 = bs; 88649bd79ccSDebojyoti Ghosh idiag += bs2; 88749bd79ccSDebojyoti Ghosh for (i = 1; i < m; i++) { 88849bd79ccSDebojyoti Ghosh v = aa + ai[i]; 88949bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 89049bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 89149bd79ccSDebojyoti Ghosh 89249bd79ccSDebojyoti Ghosh if (T) { /* b - T (Arow * x) */ 8939566063dSJacob Faibussowitsch PetscCall(PetscMemzero(w, bs * sizeof(PetscScalar))); 89449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 89549bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k]; 89649bd79ccSDebojyoti Ghosh } 89749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs, w, T, &t[i2]); 89849bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] += b[i2 + k]; 89949bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9009566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar))); 90149bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 90249bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] -= v[j] * x[vi[j] * bs + k]; 90349bd79ccSDebojyoti Ghosh } 90449bd79ccSDebojyoti Ghosh } else { 9059566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar))); 90649bd79ccSDebojyoti Ghosh } 90749bd79ccSDebojyoti Ghosh 90849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, t + i2, idiag, y); 90949bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = omega * y[j]; 91049bd79ccSDebojyoti Ghosh 91149bd79ccSDebojyoti Ghosh idiag += bs2; 91249bd79ccSDebojyoti Ghosh i2 += bs; 91349bd79ccSDebojyoti Ghosh } 91449bd79ccSDebojyoti Ghosh /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 9159566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * a->nz)); 91649bd79ccSDebojyoti Ghosh xb = t; 91749bd79ccSDebojyoti Ghosh } else xb = b; 91849bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 91949bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1); 92049bd79ccSDebojyoti Ghosh i2 = bs * (m - 1); 9219566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 92249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, x + i2); 92349bd79ccSDebojyoti Ghosh i2 -= bs; 92449bd79ccSDebojyoti Ghosh idiag -= bs2; 92549bd79ccSDebojyoti Ghosh for (i = m - 2; i >= 0; i--) { 92649bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 92749bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 92849bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 92949bd79ccSDebojyoti Ghosh 93049bd79ccSDebojyoti Ghosh if (T) { /* FIXME: This branch untested */ 9319566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 93249bd79ccSDebojyoti Ghosh /* copy all rows of x that are needed into contiguous space */ 93349bd79ccSDebojyoti Ghosh workt = work; 93449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9359566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 93649bd79ccSDebojyoti Ghosh workt += bs; 93749bd79ccSDebojyoti Ghosh } 93849bd79ccSDebojyoti Ghosh arrt = arr; 93949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9409566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 94149bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 94249bd79ccSDebojyoti Ghosh arrt += bs2; 94349bd79ccSDebojyoti Ghosh } 94449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 94549bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9469566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, t + i2, bs * sizeof(PetscScalar))); 94749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 94849bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k]; 94949bd79ccSDebojyoti Ghosh } 95049bd79ccSDebojyoti Ghosh } 95149bd79ccSDebojyoti Ghosh 95249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); /* RHS incorrect for omega != 1.0 */ 95349bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = (1.0 - omega) * x[i2 + j] + omega * y[j]; 95449bd79ccSDebojyoti Ghosh 95549bd79ccSDebojyoti Ghosh idiag -= bs2; 95649bd79ccSDebojyoti Ghosh i2 -= bs; 95749bd79ccSDebojyoti Ghosh } 9589566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz))); 95949bd79ccSDebojyoti Ghosh } 96049bd79ccSDebojyoti Ghosh its--; 96149bd79ccSDebojyoti Ghosh } 96249bd79ccSDebojyoti Ghosh while (its--) { /* FIXME: This branch not updated */ 96349bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 96449bd79ccSDebojyoti Ghosh i2 = 0; 96549bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 96649bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 9679566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar))); 96849bd79ccSDebojyoti Ghosh 96949bd79ccSDebojyoti Ghosh v = aa + ai[i]; 97049bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 97149bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 97249bd79ccSDebojyoti Ghosh workt = work; 97349bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9749566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 97549bd79ccSDebojyoti Ghosh workt += bs; 97649bd79ccSDebojyoti Ghosh } 97749bd79ccSDebojyoti Ghosh arrt = arr; 97849bd79ccSDebojyoti Ghosh if (T) { 97949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9809566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 98149bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 98249bd79ccSDebojyoti Ghosh arrt += bs2; 98349bd79ccSDebojyoti Ghosh } 98449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 98549bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 98649bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9879566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 98849bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 98949bd79ccSDebojyoti Ghosh arrt += bs2; 99049bd79ccSDebojyoti Ghosh } 99149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 99249bd79ccSDebojyoti Ghosh } 9939566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, w, bs * sizeof(PetscScalar))); 99449bd79ccSDebojyoti Ghosh 99549bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 99649bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 99749bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 99849bd79ccSDebojyoti Ghosh workt = work; 99949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10009566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 100149bd79ccSDebojyoti Ghosh workt += bs; 100249bd79ccSDebojyoti Ghosh } 100349bd79ccSDebojyoti Ghosh arrt = arr; 100449bd79ccSDebojyoti Ghosh if (T) { 100549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10069566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 100749bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 100849bd79ccSDebojyoti Ghosh arrt += bs2; 100949bd79ccSDebojyoti Ghosh } 101049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 101149bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 101249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10139566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 101449bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 101549bd79ccSDebojyoti Ghosh arrt += bs2; 101649bd79ccSDebojyoti Ghosh } 101749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 101849bd79ccSDebojyoti Ghosh } 101949bd79ccSDebojyoti Ghosh 102049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 102149bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 102249bd79ccSDebojyoti Ghosh 102349bd79ccSDebojyoti Ghosh idiag += bs2; 102449bd79ccSDebojyoti Ghosh i2 += bs; 102549bd79ccSDebojyoti Ghosh } 102649bd79ccSDebojyoti Ghosh xb = t; 10279371c9d4SSatish Balay } else xb = b; 102849bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 102949bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1); 103049bd79ccSDebojyoti Ghosh i2 = bs * (m - 1); 103149bd79ccSDebojyoti Ghosh if (xb == b) { 103249bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) { 10339566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar))); 103449bd79ccSDebojyoti Ghosh 103549bd79ccSDebojyoti Ghosh v = aa + ai[i]; 103649bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 103749bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 103849bd79ccSDebojyoti Ghosh workt = work; 103949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10409566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 104149bd79ccSDebojyoti Ghosh workt += bs; 104249bd79ccSDebojyoti Ghosh } 104349bd79ccSDebojyoti Ghosh arrt = arr; 104449bd79ccSDebojyoti Ghosh if (T) { 104549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10469566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 104749bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 104849bd79ccSDebojyoti Ghosh arrt += bs2; 104949bd79ccSDebojyoti Ghosh } 105049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 105149bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 105249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10539566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 105449bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 105549bd79ccSDebojyoti Ghosh arrt += bs2; 105649bd79ccSDebojyoti Ghosh } 105749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 105849bd79ccSDebojyoti Ghosh } 105949bd79ccSDebojyoti Ghosh 106049bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 106149bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 106249bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 106349bd79ccSDebojyoti Ghosh workt = work; 106449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10659566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 106649bd79ccSDebojyoti Ghosh workt += bs; 106749bd79ccSDebojyoti Ghosh } 106849bd79ccSDebojyoti Ghosh arrt = arr; 106949bd79ccSDebojyoti Ghosh if (T) { 107049bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10719566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 107249bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 107349bd79ccSDebojyoti Ghosh arrt += bs2; 107449bd79ccSDebojyoti Ghosh } 107549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 107649bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 107749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10789566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 107949bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 108049bd79ccSDebojyoti Ghosh arrt += bs2; 108149bd79ccSDebojyoti Ghosh } 108249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 108349bd79ccSDebojyoti Ghosh } 108449bd79ccSDebojyoti Ghosh 108549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 108649bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 108749bd79ccSDebojyoti Ghosh } 108849bd79ccSDebojyoti Ghosh } else { 108949bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) { 10909566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 109149bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 109249bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 109349bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 109449bd79ccSDebojyoti Ghosh workt = work; 109549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10969566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 109749bd79ccSDebojyoti Ghosh workt += bs; 109849bd79ccSDebojyoti Ghosh } 109949bd79ccSDebojyoti Ghosh arrt = arr; 110049bd79ccSDebojyoti Ghosh if (T) { 110149bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 11029566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 110349bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 110449bd79ccSDebojyoti Ghosh arrt += bs2; 110549bd79ccSDebojyoti Ghosh } 110649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 110749bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 110849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 11099566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 111049bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 111149bd79ccSDebojyoti Ghosh arrt += bs2; 111249bd79ccSDebojyoti Ghosh } 111349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 111449bd79ccSDebojyoti Ghosh } 111549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 111649bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 111749bd79ccSDebojyoti Ghosh } 111849bd79ccSDebojyoti Ghosh } 11199566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz))); 112049bd79ccSDebojyoti Ghosh } 112149bd79ccSDebojyoti Ghosh } 112249bd79ccSDebojyoti Ghosh 1123d0609cedSBarry Smith PetscCall(VecRestoreArray(xx, &x)); 11249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 11253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112649bd79ccSDebojyoti Ghosh } 112749bd79ccSDebojyoti Ghosh 112849bd79ccSDebojyoti Ghosh /*===================================================================================*/ 112949bd79ccSDebojyoti Ghosh 1130ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_MPIKAIJ(Mat A, Vec xx, Vec yy, Vec zz) 1131d71ae5a4SJacob Faibussowitsch { 113249bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 113349bd79ccSDebojyoti Ghosh 113449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 113549bd79ccSDebojyoti Ghosh if (!yy) { 11369566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 113749bd79ccSDebojyoti Ghosh } else { 11389566063dSJacob Faibussowitsch PetscCall(VecCopy(yy, zz)); 113949bd79ccSDebojyoti Ghosh } 11409566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */ 114149bd79ccSDebojyoti Ghosh /* start the scatter */ 11429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); 11439566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->multadd)(b->AIJ, xx, zz, zz)); 11449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); 11459566063dSJacob Faibussowitsch PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, zz, zz)); 11463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114749bd79ccSDebojyoti Ghosh } 114849bd79ccSDebojyoti Ghosh 1149ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_MPIKAIJ(Mat A, Vec xx, Vec yy) 1150d71ae5a4SJacob Faibussowitsch { 115149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1152f3fa974cSJacob Faibussowitsch PetscCall(MatMultAdd_MPIKAIJ(A, xx, NULL, yy)); 11533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115449bd79ccSDebojyoti Ghosh } 115549bd79ccSDebojyoti Ghosh 1156ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A, const PetscScalar **values) 1157d71ae5a4SJacob Faibussowitsch { 115849bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 115949bd79ccSDebojyoti Ghosh 116049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 11619566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ is up to date. */ 11629566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->invertblockdiagonal)(b->AIJ, values)); 11633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116449bd79ccSDebojyoti Ghosh } 116549bd79ccSDebojyoti Ghosh 1166ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values) 1167d71ae5a4SJacob Faibussowitsch { 116849bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 11693ba16761SJacob Faibussowitsch PetscBool diag = PETSC_FALSE; 117049bd79ccSDebojyoti Ghosh PetscInt nzaij, nz, *colsaij, *idx, i, j, p = b->p, q = b->q, r = row / p, s = row % p, c; 117149bd79ccSDebojyoti Ghosh PetscScalar *vaij, *v, *S = b->S, *T = b->T; 117249bd79ccSDebojyoti Ghosh 117349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 117428b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active"); 117549bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 1176aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row); 117749bd79ccSDebojyoti Ghosh 117849bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 117949bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 118049bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 118149bd79ccSDebojyoti Ghosh if (values) *values = NULL; 11823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118349bd79ccSDebojyoti Ghosh } 118449bd79ccSDebojyoti Ghosh 118549bd79ccSDebojyoti Ghosh if (T || b->isTI) { 11869566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(b->AIJ, r, &nzaij, &colsaij, &vaij)); 118749bd79ccSDebojyoti Ghosh c = nzaij; 118849bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 118949bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 119049bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 119149bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 119249bd79ccSDebojyoti Ghosh c = i; 119349bd79ccSDebojyoti Ghosh } 119449bd79ccSDebojyoti Ghosh } 119549bd79ccSDebojyoti Ghosh } else nzaij = c = 0; 119649bd79ccSDebojyoti Ghosh 119749bd79ccSDebojyoti Ghosh /* calculate size of row */ 119849bd79ccSDebojyoti Ghosh nz = 0; 119949bd79ccSDebojyoti Ghosh if (S) nz += q; 120049bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (nzaij - 1) * q : nzaij * q); 120149bd79ccSDebojyoti Ghosh 120249bd79ccSDebojyoti Ghosh if (cols || values) { 12039566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz, &idx, nz, &v)); 120438822f9dSRichard Tran Mills for (i = 0; i < q; i++) { 120538822f9dSRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 120638822f9dSRichard Tran Mills v[i] = 0.0; 120738822f9dSRichard Tran Mills } 120849bd79ccSDebojyoti Ghosh if (b->isTI) { 120949bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 121049bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 121149bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j; 121249bd79ccSDebojyoti Ghosh v[i * q + j] = (j == s ? vaij[i] : 0); 121349bd79ccSDebojyoti Ghosh } 121449bd79ccSDebojyoti Ghosh } 121549bd79ccSDebojyoti Ghosh } else if (T) { 121649bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 121749bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 121849bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j; 121949bd79ccSDebojyoti Ghosh v[i * q + j] = vaij[i] * T[s + j * p]; 122049bd79ccSDebojyoti Ghosh } 122149bd79ccSDebojyoti Ghosh } 122249bd79ccSDebojyoti Ghosh } 122349bd79ccSDebojyoti Ghosh if (S) { 122449bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 122549bd79ccSDebojyoti Ghosh idx[c * q + j] = r * q + j; 122649bd79ccSDebojyoti Ghosh v[c * q + j] += S[s + j * p]; 122749bd79ccSDebojyoti Ghosh } 122849bd79ccSDebojyoti Ghosh } 122949bd79ccSDebojyoti Ghosh } 123049bd79ccSDebojyoti Ghosh 123149bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 123249bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 123349bd79ccSDebojyoti Ghosh if (values) *values = v; 12343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123549bd79ccSDebojyoti Ghosh } 123649bd79ccSDebojyoti Ghosh 1237ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 1238d71ae5a4SJacob Faibussowitsch { 123949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1240cb4a9cd9SHong Zhang if (nz) *nz = 0; 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 136211a5261eSBarry Smith MatCreateKAIJ - Creates a matrix of type `MATKAIJ` to be used for matrices of the following form: 136349bd79ccSDebojyoti Ghosh 136449bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 136549bd79ccSDebojyoti Ghosh 136649bd79ccSDebojyoti Ghosh where 136749bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 136849bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 136911a5261eSBarry Smith A is an `MATAIJ` (n \times n) matrix 137049bd79ccSDebojyoti Ghosh I is the identity matrix 137149bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 137249bd79ccSDebojyoti Ghosh 1373*2ef1f0ffSBarry Smith `S` and `T` are always stored independently on all processes as `PetscScalar` arrays in column-major format. 137449bd79ccSDebojyoti Ghosh 137549bd79ccSDebojyoti Ghosh Collective 137649bd79ccSDebojyoti Ghosh 137749bd79ccSDebojyoti Ghosh Input Parameters: 137811a5261eSBarry Smith + A - the `MATAIJ` matrix 1379*2ef1f0ffSBarry Smith . p - number of rows in `S` and `T` 1380*2ef1f0ffSBarry Smith . q - number of columns in `S` and `T` 1381*2ef1f0ffSBarry Smith . S - the `S` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major) 1382*2ef1f0ffSBarry Smith - T - the `T` matrix (can be `NULL`), stored as a `PetscScalar` array (column-major) 138349bd79ccSDebojyoti Ghosh 138449bd79ccSDebojyoti Ghosh Output Parameter: 138511a5261eSBarry Smith . kaij - the new `MATKAIJ` matrix 138649bd79ccSDebojyoti Ghosh 1387*2ef1f0ffSBarry Smith Level: advanced 1388*2ef1f0ffSBarry Smith 1389d3b90ce1SRichard Tran Mills Notes: 139011a5261eSBarry 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. 139149bd79ccSDebojyoti Ghosh 139211a5261eSBarry Smith Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix. 139311a5261eSBarry Smith 139411a5261eSBarry Smith Developer Note: 139511a5261eSBarry Smith In the `MATMPIKAIJ` case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state 1396*2ef1f0ffSBarry Smith of the AIJ matrix 'A' that describes the blockwise action of the `MATMPIKAIJ` matrix and, if the object state has changed, lazily 1397*2ef1f0ffSBarry Smith rebuilding 'AIJ' and 'OAIJ' just before executing operations with the `MATMPIKAIJ` matrix. If new types of operations are added, 1398761d359dSRichard Tran Mills routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine). 1399761d359dSRichard Tran Mills 1400*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MATKAIJ` 140149bd79ccSDebojyoti Ghosh @*/ 1402d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateKAIJ(Mat A, PetscInt p, PetscInt q, const PetscScalar S[], const PetscScalar T[], Mat *kaij) 1403d71ae5a4SJacob Faibussowitsch { 140449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 14059566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), kaij)); 14069566063dSJacob Faibussowitsch PetscCall(MatSetType(*kaij, MATKAIJ)); 14079566063dSJacob Faibussowitsch PetscCall(MatKAIJSetAIJ(*kaij, A)); 14089566063dSJacob Faibussowitsch PetscCall(MatKAIJSetS(*kaij, p, q, S)); 14099566063dSJacob Faibussowitsch PetscCall(MatKAIJSetT(*kaij, p, q, T)); 14109566063dSJacob Faibussowitsch PetscCall(MatSetUp(*kaij)); 14113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14120567c835SRichard Tran Mills } 141349bd79ccSDebojyoti Ghosh 14140567c835SRichard Tran Mills /*MC 14155881e567SRichard Tran Mills MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form 14165881e567SRichard Tran Mills [I \otimes S + A \otimes T], 14170567c835SRichard Tran Mills where 1418*2ef1f0ffSBarry Smith .vb 14195881e567SRichard Tran Mills S is a dense (p \times q) matrix, 14205881e567SRichard Tran Mills T is a dense (p \times q) matrix, 14215881e567SRichard Tran Mills A is an AIJ (n \times n) matrix, 14225881e567SRichard Tran Mills and I is the identity matrix. 1423*2ef1f0ffSBarry Smith .ve 14245881e567SRichard Tran Mills The resulting matrix is (np \times nq). 14250567c835SRichard Tran Mills 142611a5261eSBarry Smith S and T are always stored independently on all processes as `PetscScalar` arrays in column-major format. 14270567c835SRichard Tran Mills 1428*2ef1f0ffSBarry Smith Level: advanced 1429*2ef1f0ffSBarry Smith 143011a5261eSBarry Smith Note: 14315881e567SRichard 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, 14325881e567SRichard Tran Mills where x and b are column vectors containing the row-major representations of X and B. 14335881e567SRichard Tran Mills 1434*2ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MatCreateKAIJ()` 14350567c835SRichard Tran Mills M*/ 14360567c835SRichard Tran Mills 1437d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A) 1438d71ae5a4SJacob Faibussowitsch { 14390567c835SRichard Tran Mills Mat_MPIKAIJ *b; 14400567c835SRichard Tran Mills PetscMPIInt size; 14410567c835SRichard Tran Mills 14420567c835SRichard Tran Mills PetscFunctionBegin; 14434dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 14440567c835SRichard Tran Mills A->data = (void *)b; 14450567c835SRichard Tran Mills 14469566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps))); 14470567c835SRichard Tran Mills 1448f4259b30SLisandro Dalcin b->w = NULL; 14499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 14500567c835SRichard Tran Mills if (size == 1) { 14519566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQKAIJ)); 14520567c835SRichard Tran Mills A->ops->destroy = MatDestroy_SeqKAIJ; 1453bb6fb833SRichard Tran Mills A->ops->mult = MatMult_SeqKAIJ; 1454bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_SeqKAIJ; 1455bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ; 14560567c835SRichard Tran Mills A->ops->getrow = MatGetRow_SeqKAIJ; 14570567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_SeqKAIJ; 14580567c835SRichard Tran Mills A->ops->sor = MatSOR_SeqKAIJ; 14599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", MatConvert_KAIJ_AIJ)); 14600567c835SRichard Tran Mills } else { 14619566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMPIKAIJ)); 14620567c835SRichard Tran Mills A->ops->destroy = MatDestroy_MPIKAIJ; 1463bb6fb833SRichard Tran Mills A->ops->mult = MatMult_MPIKAIJ; 1464bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_MPIKAIJ; 1465bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ; 14660567c835SRichard Tran Mills A->ops->getrow = MatGetRow_MPIKAIJ; 14670567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_MPIKAIJ; 14689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", MatGetDiagonalBlock_MPIKAIJ)); 14699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", MatConvert_KAIJ_AIJ)); 14700567c835SRichard Tran Mills } 147106d61a37SPierre Jolivet A->ops->setup = MatSetUp_KAIJ; 147206d61a37SPierre Jolivet A->ops->view = MatView_KAIJ; 14730567c835SRichard Tran Mills A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ; 14743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147549bd79ccSDebojyoti Ghosh } 1476