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 4611a5261eSBarry Smith .seealso: `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 6811a5261eSBarry Smith MatKAIJGetS - Get the S matrix describing the shift action of the `MATKAIJ` matrix 6949bd79ccSDebojyoti Ghosh 700567c835SRichard Tran Mills Not Collective; the entire S is stored and returned independently on all processes. 7149bd79ccSDebojyoti Ghosh 7249bd79ccSDebojyoti Ghosh Input Parameter: 7311a5261eSBarry Smith . A - the `MATKAIJ` matrix 7449bd79ccSDebojyoti Ghosh 75a5b5c723SRichard Tran Mills Output Parameters: 76a5b5c723SRichard Tran Mills + m - the number of rows in S 77a5b5c723SRichard Tran Mills . n - the number of columns in S 78a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 7949bd79ccSDebojyoti Ghosh 8011a5261eSBarry Smith Note: 8111a5261eSBarry Smith All output parameters are optional (pass NULL if not desired) 8231a97b9aSRichard Tran Mills 8349bd79ccSDebojyoti Ghosh Level: advanced 8449bd79ccSDebojyoti Ghosh 8511a5261eSBarry Smith .seealso: `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 9811a5261eSBarry Smith MatKAIJGetSRead - Get a read-only pointer to the S matrix describing the shift action of the `MATKAIJ` matrix 99a5b5c723SRichard Tran Mills 100a5b5c723SRichard Tran Mills 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: 106a5b5c723SRichard Tran Mills + m - the number of rows in S 107a5b5c723SRichard Tran Mills . 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 11011a5261eSBarry Smith Note: 11111a5261eSBarry Smith All output parameters are optional (pass NULL if not desired) 112a5b5c723SRichard Tran Mills 113a5b5c723SRichard Tran Mills Level: advanced 114a5b5c723SRichard Tran Mills 11511a5261eSBarry Smith .seealso: `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 130a5b5c723SRichard Tran Mills Not collective 131a5b5c723SRichard Tran Mills 132a5b5c723SRichard Tran Mills Input Parameter: 13311a5261eSBarry Smith . A - the `MATKAIJ` matrix 134a5b5c723SRichard Tran Mills 135a5b5c723SRichard Tran Mills Output Parameter: 13611a5261eSBarry Smith . S - location of pointer to array obtained with `MatKAIJGetS()` 137a5b5c723SRichard Tran Mills 13811a5261eSBarry Smith Note: 13911a5261eSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 140a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 141a5b5c723SRichard Tran Mills 142a5b5c723SRichard Tran Mills Level: advanced 14311a5261eSBarry Smith 14411a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`, `MatKAIJRestoreSRead()` 145a5b5c723SRichard Tran Mills @*/ 146d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreS(Mat A, PetscScalar **S) 147d71ae5a4SJacob Faibussowitsch { 148a5b5c723SRichard Tran Mills PetscFunctionBegin; 149a5b5c723SRichard Tran Mills if (S) *S = NULL; 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)A)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152a5b5c723SRichard Tran Mills } 153a5b5c723SRichard Tran Mills 154a5b5c723SRichard Tran Mills /*@C 15511a5261eSBarry Smith MatKAIJRestoreSRead - Restore array obtained with `MatKAIJGetSRead()` 156a5b5c723SRichard Tran Mills 157a5b5c723SRichard Tran Mills Not collective 158a5b5c723SRichard Tran Mills 159a5b5c723SRichard Tran Mills Input Parameter: 16011a5261eSBarry Smith . A - the `MATKAIJ` matrix 161a5b5c723SRichard Tran Mills 162a5b5c723SRichard Tran Mills Output Parameter: 16311a5261eSBarry Smith . S - location of pointer to array obtained with `MatKAIJGetS()` 164a5b5c723SRichard Tran Mills 16511a5261eSBarry Smith Note: 16611a5261eSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 167a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 168a5b5c723SRichard Tran Mills 169a5b5c723SRichard Tran Mills Level: advanced 17011a5261eSBarry Smith 17111a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetSRead()`, `MatKAIJRestoreSRead()` 172a5b5c723SRichard Tran Mills @*/ 173d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreSRead(Mat A, const PetscScalar **S) 174d71ae5a4SJacob Faibussowitsch { 175a5b5c723SRichard Tran Mills PetscFunctionBegin; 176a5b5c723SRichard Tran Mills if (S) *S = NULL; 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17849bd79ccSDebojyoti Ghosh } 17949bd79ccSDebojyoti Ghosh 18049bd79ccSDebojyoti Ghosh /*@C 18111a5261eSBarry Smith MatKAIJGetT - Get the transformation matrix T associated with the `MATKAIJ` matrix 18249bd79ccSDebojyoti Ghosh 1830567c835SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 18449bd79ccSDebojyoti Ghosh 18549bd79ccSDebojyoti Ghosh Input Parameter: 18611a5261eSBarry Smith . A - the `MATKAIJ` matrix 18749bd79ccSDebojyoti Ghosh 188d8d19677SJose E. Roman Output Parameters: 189a5b5c723SRichard Tran Mills + m - the number of rows in T 190a5b5c723SRichard Tran Mills . n - the number of columns in T 191a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 19249bd79ccSDebojyoti Ghosh 19311a5261eSBarry Smith Note: 19411a5261eSBarry Smith All output parameters are optional (pass NULL if not desired) 19531a97b9aSRichard Tran Mills 19649bd79ccSDebojyoti Ghosh Level: advanced 19749bd79ccSDebojyoti Ghosh 19811a5261eSBarry Smith .seealso: `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 19949bd79ccSDebojyoti Ghosh @*/ 200d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetT(Mat A, PetscInt *m, PetscInt *n, PetscScalar **T) 201d71ae5a4SJacob Faibussowitsch { 20249bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 20349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 204a5b5c723SRichard Tran Mills if (m) *m = b->p; 205a5b5c723SRichard Tran Mills if (n) *n = b->q; 206a5b5c723SRichard Tran Mills if (T) *T = b->T; 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 208a5b5c723SRichard Tran Mills } 209a5b5c723SRichard Tran Mills 210a5b5c723SRichard Tran Mills /*@C 21111a5261eSBarry Smith MatKAIJGetTRead - Get a read-only pointer to the transformation matrix T associated with the `MATKAIJ` matrix 212a5b5c723SRichard Tran Mills 213a5b5c723SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 214a5b5c723SRichard Tran Mills 215a5b5c723SRichard Tran Mills Input Parameter: 21611a5261eSBarry Smith . A - the `MATKAIJ` matrix 217a5b5c723SRichard Tran Mills 218d8d19677SJose E. Roman Output Parameters: 219a5b5c723SRichard Tran Mills + m - the number of rows in T 220a5b5c723SRichard Tran Mills . n - the number of columns in T 221a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 222a5b5c723SRichard Tran Mills 22311a5261eSBarry Smith Note: All output parameters are optional (pass NULL if not desired) 224a5b5c723SRichard Tran Mills 225a5b5c723SRichard Tran Mills Level: advanced 226a5b5c723SRichard Tran Mills 22711a5261eSBarry Smith .seealso: `MATKAIJ`, `MatCreateKAIJ()`, `MatGetBlockSizes()` 228a5b5c723SRichard Tran Mills @*/ 229d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetTRead(Mat A, PetscInt *m, PetscInt *n, const PetscScalar **T) 230d71ae5a4SJacob Faibussowitsch { 231a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 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 242a5b5c723SRichard Tran Mills Not collective 243a5b5c723SRichard Tran Mills 244a5b5c723SRichard Tran Mills Input Parameter: 24511a5261eSBarry Smith . A - the `MATKAIJ` matrix 246a5b5c723SRichard Tran Mills 247a5b5c723SRichard Tran Mills Output Parameter: 24811a5261eSBarry Smith . T - location of pointer to array obtained with `MatKAIJGetS()` 249a5b5c723SRichard Tran Mills 25011a5261eSBarry Smith Note: 25111a5261eSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 252a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 253a5b5c723SRichard Tran Mills 254a5b5c723SRichard Tran Mills Level: advanced 25511a5261eSBarry Smith 25611a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`, `MatKAIJRestoreTRead()` 257a5b5c723SRichard Tran Mills @*/ 258d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreT(Mat A, PetscScalar **T) 259d71ae5a4SJacob Faibussowitsch { 260a5b5c723SRichard Tran Mills PetscFunctionBegin; 261a5b5c723SRichard Tran Mills if (T) *T = NULL; 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectStateIncrease((PetscObject)A)); 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 264a5b5c723SRichard Tran Mills } 265a5b5c723SRichard Tran Mills 266a5b5c723SRichard Tran Mills /*@C 26711a5261eSBarry Smith MatKAIJRestoreTRead - Restore array obtained with `MatKAIJGetTRead()` 268a5b5c723SRichard Tran Mills 269a5b5c723SRichard Tran Mills Not collective 270a5b5c723SRichard Tran Mills 271a5b5c723SRichard Tran Mills Input Parameter: 27211a5261eSBarry Smith . A - the `MATKAIJ` matrix 273a5b5c723SRichard Tran Mills 274a5b5c723SRichard Tran Mills Output Parameter: 27511a5261eSBarry Smith . T - location of pointer to array obtained with `MatKAIJGetS()` 276a5b5c723SRichard Tran Mills 27711a5261eSBarry Smith Note: 27811a5261eSBarry Smith This routine zeros the array pointer to prevent accidental reuse after it has been restored. 279a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 280a5b5c723SRichard Tran Mills 281a5b5c723SRichard Tran Mills Level: advanced 28211a5261eSBarry Smith 28311a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJGetTRead()`, `MatKAIJRestoreTRead()` 284a5b5c723SRichard Tran Mills @*/ 285d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJRestoreTRead(Mat A, const PetscScalar **T) 286d71ae5a4SJacob Faibussowitsch { 287a5b5c723SRichard Tran Mills PetscFunctionBegin; 288a5b5c723SRichard Tran Mills if (T) *T = NULL; 2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29049bd79ccSDebojyoti Ghosh } 29149bd79ccSDebojyoti Ghosh 2920567c835SRichard Tran Mills /*@ 29311a5261eSBarry Smith MatKAIJSetAIJ - Set the `MATAIJ` matrix describing the blockwise action of the `MATKAIJ` matrix 2940567c835SRichard Tran Mills 29511a5261eSBarry Smith Logically Collective; if the `MATAIJ` matrix is parallel, the `MATKAIJ` matrix is also parallel 2960567c835SRichard Tran Mills 2970567c835SRichard Tran Mills Input Parameters: 29811a5261eSBarry Smith + A - the `MATKAIJ` matrix 29911a5261eSBarry Smith - B - the `MATAIJ` matrix 3000567c835SRichard Tran Mills 30115b9d025SRichard Tran Mills Notes: 30211a5261eSBarry 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. 30311a5261eSBarry Smith 30411a5261eSBarry Smith Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix. 30515b9d025SRichard Tran Mills 3060567c835SRichard Tran Mills Level: advanced 3070567c835SRichard Tran Mills 30811a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()` 3090567c835SRichard Tran Mills @*/ 310d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetAIJ(Mat A, Mat B) 311d71ae5a4SJacob Faibussowitsch { 3120567c835SRichard Tran Mills PetscMPIInt size; 31306d61a37SPierre Jolivet PetscBool flg; 3140567c835SRichard Tran Mills 3150567c835SRichard Tran Mills PetscFunctionBegin; 3169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 3170567c835SRichard Tran Mills if (size == 1) { 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &flg)); 31928b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)B), PETSC_ERR_SUP, "MatKAIJSetAIJ() with MATSEQKAIJ does not support %s as the AIJ mat", ((PetscObject)B)->type_name); 3200567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 3210567c835SRichard Tran Mills a->AIJ = B; 3220567c835SRichard Tran Mills } else { 3230567c835SRichard Tran Mills Mat_MPIKAIJ *a = (Mat_MPIKAIJ *)A->data; 3240567c835SRichard Tran Mills a->A = B; 3250567c835SRichard Tran Mills } 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)B)); 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3280567c835SRichard Tran Mills } 3290567c835SRichard Tran Mills 3300567c835SRichard Tran Mills /*@C 33111a5261eSBarry Smith MatKAIJSetS - Set the S matrix describing the shift action of the `MATKAIJ` matrix 3320567c835SRichard Tran Mills 3330567c835SRichard Tran Mills Logically Collective; the entire S is stored independently on all processes. 3340567c835SRichard Tran Mills 3350567c835SRichard Tran Mills Input Parameters: 33611a5261eSBarry Smith + A - the `MATKAIJ` matrix 3370567c835SRichard Tran Mills . p - the number of rows in S 3380567c835SRichard Tran Mills . q - the number of columns in S 3390567c835SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 3400567c835SRichard Tran Mills 34111a5261eSBarry Smith Notes: 34211a5261eSBarry Smith The dimensions p and q must match those of the transformation matrix T associated with the `MATKAIJ` matrix. 34311a5261eSBarry Smith 34488f48298SRichard Tran Mills The S matrix is copied, so the user can destroy this array. 3450567c835SRichard Tran Mills 34611a5261eSBarry Smith Level: advanced 3470567c835SRichard Tran Mills 34811a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJSetT()`, `MatKAIJSetAIJ()` 3490567c835SRichard Tran Mills @*/ 350d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetS(Mat A, PetscInt p, PetscInt q, const PetscScalar S[]) 351d71ae5a4SJacob Faibussowitsch { 3520567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 3530567c835SRichard Tran Mills 3540567c835SRichard Tran Mills PetscFunctionBegin; 3559566063dSJacob Faibussowitsch PetscCall(PetscFree(a->S)); 3560567c835SRichard Tran Mills if (S) { 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p * q, &a->S)); 3589566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(a->S, S, p * q * sizeof(PetscScalar))); 3590567c835SRichard Tran Mills } else a->S = NULL; 3600567c835SRichard Tran Mills 3610567c835SRichard Tran Mills a->p = p; 3620567c835SRichard Tran Mills a->q = q; 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3640567c835SRichard Tran Mills } 3650567c835SRichard Tran Mills 3660567c835SRichard Tran Mills /*@C 367910cf402Sprj- MatKAIJGetScaledIdentity - Check if both S and T are scaled identities. 368910cf402Sprj- 369910cf402Sprj- Logically Collective. 370910cf402Sprj- 371910cf402Sprj- Input Parameter: 37211a5261eSBarry Smith . A - the `MATKAIJ` matrix 373910cf402Sprj- 374910cf402Sprj- Output Parameter: 375910cf402Sprj- . identity - the Boolean value 376910cf402Sprj- 377910cf402Sprj- Level: Advanced 378910cf402Sprj- 37911a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetS()`, `MatKAIJGetT()` 380910cf402Sprj- @*/ 381d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJGetScaledIdentity(Mat A, PetscBool *identity) 382d71ae5a4SJacob Faibussowitsch { 383910cf402Sprj- Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 384910cf402Sprj- PetscInt i, j; 385910cf402Sprj- 386910cf402Sprj- PetscFunctionBegin; 387910cf402Sprj- if (a->p != a->q) { 388910cf402Sprj- *identity = PETSC_FALSE; 3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 390910cf402Sprj- } else *identity = PETSC_TRUE; 391910cf402Sprj- if (!a->isTI || a->S) { 392910cf402Sprj- for (i = 0; i < a->p && *identity; i++) { 393910cf402Sprj- for (j = 0; j < a->p && *identity; j++) { 394910cf402Sprj- if (i != j) { 395910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 396910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i + j * a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 397910cf402Sprj- } else { 398910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i * (a->p + 1)] - a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 399910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i * (a->p + 1)] - a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 400910cf402Sprj- } 401910cf402Sprj- } 402910cf402Sprj- } 403910cf402Sprj- } 4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 405910cf402Sprj- } 406910cf402Sprj- 407910cf402Sprj- /*@C 40811a5261eSBarry Smith MatKAIJSetT - Set the transformation matrix T associated with the `MATKAIJ` matrix 4090567c835SRichard Tran Mills 4100567c835SRichard Tran Mills Logically Collective; the entire T is stored independently on all processes. 4110567c835SRichard Tran Mills 4120567c835SRichard Tran Mills Input Parameters: 41311a5261eSBarry Smith + A - the `MATKAIJ` matrix 4140567c835SRichard Tran Mills . p - the number of rows in S 4150567c835SRichard Tran Mills . q - the number of columns in S 4160567c835SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 4170567c835SRichard Tran Mills 41811a5261eSBarry Smith Notes: 41911a5261eSBarry Smith The dimensions p and q must match those of the shift matrix S associated with the `MATKAIJ` matrix. 42011a5261eSBarry Smith 42188f48298SRichard Tran Mills The T matrix is copied, so the user can destroy this array. 4220567c835SRichard Tran Mills 4230567c835SRichard Tran Mills Level: Advanced 4240567c835SRichard Tran Mills 42511a5261eSBarry Smith .seealso: `MATKAIJ`, `MatKAIJGetT()`, `MatKAIJSetS()`, `MatKAIJSetAIJ()` 4260567c835SRichard Tran Mills @*/ 427d71ae5a4SJacob Faibussowitsch PetscErrorCode MatKAIJSetT(Mat A, PetscInt p, PetscInt q, const PetscScalar T[]) 428d71ae5a4SJacob Faibussowitsch { 4290567c835SRichard Tran Mills PetscInt i, j; 4300567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 4310567c835SRichard Tran Mills PetscBool isTI = PETSC_FALSE; 4320567c835SRichard Tran Mills 4330567c835SRichard Tran Mills PetscFunctionBegin; 4340567c835SRichard Tran Mills /* check if T is an identity matrix */ 4350567c835SRichard Tran Mills if (T && (p == q)) { 4360567c835SRichard Tran Mills isTI = PETSC_TRUE; 4370567c835SRichard Tran Mills for (i = 0; i < p; i++) { 4380567c835SRichard Tran Mills for (j = 0; j < q; j++) { 4390567c835SRichard Tran Mills if (i == j) { 4400567c835SRichard Tran Mills /* diagonal term must be 1 */ 4410567c835SRichard Tran Mills if (T[i + j * p] != 1.0) isTI = PETSC_FALSE; 4420567c835SRichard Tran Mills } else { 4430567c835SRichard Tran Mills /* off-diagonal term must be 0 */ 4440567c835SRichard Tran Mills if (T[i + j * p] != 0.0) isTI = PETSC_FALSE; 4450567c835SRichard Tran Mills } 4460567c835SRichard Tran Mills } 4470567c835SRichard Tran Mills } 4480567c835SRichard Tran Mills } 4490567c835SRichard Tran Mills a->isTI = isTI; 4500567c835SRichard Tran Mills 4519566063dSJacob Faibussowitsch PetscCall(PetscFree(a->T)); 4520567c835SRichard Tran Mills if (T && (!isTI)) { 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p * q, &a->T)); 4549566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(a->T, T, p * q * sizeof(PetscScalar))); 45550d19d74SRichard Tran Mills } else a->T = NULL; 4560567c835SRichard Tran Mills 4570567c835SRichard Tran Mills a->p = p; 4580567c835SRichard Tran Mills a->q = q; 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4600567c835SRichard Tran Mills } 4610567c835SRichard Tran Mills 462*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_SeqKAIJ(Mat A) 463d71ae5a4SJacob Faibussowitsch { 46449bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 46549bd79ccSDebojyoti Ghosh 46649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->AIJ)); 4689566063dSJacob Faibussowitsch PetscCall(PetscFree(b->S)); 4699566063dSJacob Faibussowitsch PetscCall(PetscFree(b->T)); 4709566063dSJacob Faibussowitsch PetscCall(PetscFree(b->ibdiag)); 4719566063dSJacob Faibussowitsch PetscCall(PetscFree5(b->sor.w, b->sor.y, b->sor.work, b->sor.t, b->sor.arr)); 4729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", NULL)); 4739566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47549bd79ccSDebojyoti Ghosh } 47649bd79ccSDebojyoti Ghosh 477*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatKAIJ_build_AIJ_OAIJ(Mat A) 478d71ae5a4SJacob Faibussowitsch { 479e0e5a793SRichard Tran Mills Mat_MPIKAIJ *a; 480e0e5a793SRichard Tran Mills Mat_MPIAIJ *mpiaij; 481e0e5a793SRichard Tran Mills PetscScalar *T; 482e0e5a793SRichard Tran Mills PetscInt i, j; 483e0e5a793SRichard Tran Mills PetscObjectState state; 484e0e5a793SRichard Tran Mills 485e0e5a793SRichard Tran Mills PetscFunctionBegin; 486e0e5a793SRichard Tran Mills a = (Mat_MPIKAIJ *)A->data; 487e0e5a793SRichard Tran Mills mpiaij = (Mat_MPIAIJ *)a->A->data; 488e0e5a793SRichard Tran Mills 4899566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)a->A, &state)); 490e0e5a793SRichard Tran Mills if (state == a->state) { 491e0e5a793SRichard Tran Mills /* The existing AIJ and KAIJ members are up-to-date, so simply exit. */ 4923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 493e0e5a793SRichard Tran Mills } else { 4949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&a->AIJ)); 4959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&a->OAIJ)); 496e0e5a793SRichard Tran Mills if (a->isTI) { 497e0e5a793SRichard Tran Mills /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL. 498e0e5a793SRichard Tran Mills * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will 499e0e5a793SRichard Tran Mills * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix 500e0e5a793SRichard Tran Mills * to pass in. */ 5019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->p * a->q, &T)); 502e0e5a793SRichard Tran Mills for (i = 0; i < a->p; i++) { 503e0e5a793SRichard Tran Mills for (j = 0; j < a->q; j++) { 504e0e5a793SRichard Tran Mills if (i == j) T[i + j * a->p] = 1.0; 505e0e5a793SRichard Tran Mills else T[i + j * a->p] = 0.0; 506e0e5a793SRichard Tran Mills } 507e0e5a793SRichard Tran Mills } 508e0e5a793SRichard Tran Mills } else T = a->T; 5099566063dSJacob Faibussowitsch PetscCall(MatCreateKAIJ(mpiaij->A, a->p, a->q, a->S, T, &a->AIJ)); 5109566063dSJacob Faibussowitsch PetscCall(MatCreateKAIJ(mpiaij->B, a->p, a->q, NULL, T, &a->OAIJ)); 5111baa6e33SBarry Smith if (a->isTI) PetscCall(PetscFree(T)); 512e0e5a793SRichard Tran Mills a->state = state; 513e0e5a793SRichard Tran Mills } 514e0e5a793SRichard Tran Mills 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516e0e5a793SRichard Tran Mills } 517e0e5a793SRichard Tran Mills 518*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSetUp_KAIJ(Mat A) 519d71ae5a4SJacob Faibussowitsch { 5200567c835SRichard Tran Mills PetscInt n; 5210567c835SRichard Tran Mills PetscMPIInt size; 5220567c835SRichard Tran Mills Mat_SeqKAIJ *seqkaij = (Mat_SeqKAIJ *)A->data; 5230567c835SRichard Tran Mills 52449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 5260567c835SRichard Tran Mills if (size == 1) { 5279566063dSJacob 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)); 5289566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p)); 5299566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q)); 5309566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5319566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5320567c835SRichard Tran Mills } else { 5330567c835SRichard Tran Mills Mat_MPIKAIJ *a; 5340567c835SRichard Tran Mills Mat_MPIAIJ *mpiaij; 5350567c835SRichard Tran Mills IS from, to; 5360567c835SRichard Tran Mills Vec gvec; 5370567c835SRichard Tran Mills 5380567c835SRichard Tran Mills a = (Mat_MPIKAIJ *)A->data; 539d3f912faSRichard Tran Mills mpiaij = (Mat_MPIAIJ *)a->A->data; 5409566063dSJacob 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)); 5419566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->rmap, seqkaij->p)); 5429566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(A->cmap, seqkaij->q)); 5439566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5449566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5450567c835SRichard Tran Mills 5469566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); 5470567c835SRichard Tran Mills 5489566063dSJacob Faibussowitsch PetscCall(VecGetSize(mpiaij->lvec, &n)); 5499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &a->w)); 5509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(a->w, n * a->q, n * a->q)); 5519566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(a->w, a->q)); 5529566063dSJacob Faibussowitsch PetscCall(VecSetType(a->w, VECSEQ)); 5530567c835SRichard Tran Mills 5540567c835SRichard Tran Mills /* create two temporary Index sets for build scatter gather */ 5559566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)a->A), a->q, n, mpiaij->garray, PETSC_COPY_VALUES, &from)); 5569566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, n * a->q, 0, 1, &to)); 5570567c835SRichard Tran Mills 5580567c835SRichard Tran Mills /* create temporary global vector to generate scatter context */ 5599566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A), a->q, a->q * a->A->cmap->n, a->q * a->A->cmap->N, NULL, &gvec)); 5600567c835SRichard Tran Mills 5610567c835SRichard Tran Mills /* generate the scatter context */ 5629566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(gvec, from, a->w, to, &a->ctx)); 5630567c835SRichard Tran Mills 5649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 5659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 5669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gvec)); 5670567c835SRichard Tran Mills } 5680567c835SRichard Tran Mills 5690567c835SRichard Tran Mills A->assembled = PETSC_TRUE; 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57149bd79ccSDebojyoti Ghosh } 57249bd79ccSDebojyoti Ghosh 573*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatView_KAIJ(Mat A, PetscViewer viewer) 574d71ae5a4SJacob Faibussowitsch { 575e6985dafSRichard Tran Mills PetscViewerFormat format; 576e6985dafSRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 57749bd79ccSDebojyoti Ghosh Mat B; 578e6985dafSRichard Tran Mills PetscInt i; 579e6985dafSRichard Tran Mills PetscBool ismpikaij; 58049bd79ccSDebojyoti Ghosh 58149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij)); 5839566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 584e6985dafSRichard Tran Mills if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) { 5859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "S and T have %" PetscInt_FMT " rows and %" PetscInt_FMT " columns\n", a->p, a->q)); 586e6985dafSRichard Tran Mills 587e6985dafSRichard Tran Mills /* Print appropriate details for S. */ 588e6985dafSRichard Tran Mills if (!a->S) { 5899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "S is NULL\n")); 590e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 5919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of S are ")); 592e6985dafSRichard Tran Mills for (i = 0; i < (a->p * a->q); i++) { 593e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 5949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->S[i]), (double)PetscImaginaryPart(a->S[i]))); 595e6985dafSRichard Tran Mills #else 5969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->S[i]))); 597e6985dafSRichard Tran Mills #endif 598e6985dafSRichard Tran Mills } 5999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 60049bd79ccSDebojyoti Ghosh } 60149bd79ccSDebojyoti Ghosh 602e6985dafSRichard Tran Mills /* Print appropriate details for T. */ 603e6985dafSRichard Tran Mills if (a->isTI) { 6049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "T is the identity matrix\n")); 605e6985dafSRichard Tran Mills } else if (!a->T) { 6069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "T is NULL\n")); 607e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 6089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Entries of T are ")); 609e6985dafSRichard Tran Mills for (i = 0; i < (a->p * a->q); i++) { 610e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 6119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e ", (double)PetscRealPart(a->T[i]), (double)PetscImaginaryPart(a->T[i]))); 612e6985dafSRichard Tran Mills #else 6139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e ", (double)PetscRealPart(a->T[i]))); 614e6985dafSRichard Tran Mills #endif 615e6985dafSRichard Tran Mills } 6169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 617e6985dafSRichard Tran Mills } 61849bd79ccSDebojyoti Ghosh 619e6985dafSRichard Tran Mills /* Now print details for the AIJ matrix, using the AIJ viewer. */ 6209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Now viewing the associated AIJ matrix:\n")); 621e6985dafSRichard Tran Mills if (ismpikaij) { 622e6985dafSRichard Tran Mills Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 6239566063dSJacob Faibussowitsch PetscCall(MatView(b->A, viewer)); 624e6985dafSRichard Tran Mills } else { 6259566063dSJacob Faibussowitsch PetscCall(MatView(a->AIJ, viewer)); 626e6985dafSRichard Tran Mills } 627e6985dafSRichard Tran Mills 628e6985dafSRichard Tran Mills } else { 629e6985dafSRichard Tran Mills /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */ 6309566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B)); 6319566063dSJacob Faibussowitsch PetscCall(MatView(B, viewer)); 6329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 633e6985dafSRichard Tran Mills } 6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63549bd79ccSDebojyoti Ghosh } 63649bd79ccSDebojyoti Ghosh 637*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatDestroy_MPIKAIJ(Mat A) 638d71ae5a4SJacob Faibussowitsch { 63949bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 64049bd79ccSDebojyoti Ghosh 64149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 6429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->AIJ)); 6439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->OAIJ)); 6449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&b->A)); 6459566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&b->ctx)); 6469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->w)); 6479566063dSJacob Faibussowitsch PetscCall(PetscFree(b->S)); 6489566063dSJacob Faibussowitsch PetscCall(PetscFree(b->T)); 6499566063dSJacob Faibussowitsch PetscCall(PetscFree(b->ibdiag)); 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", NULL)); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", NULL)); 6529566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data)); 6533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65449bd79ccSDebojyoti Ghosh } 65549bd79ccSDebojyoti Ghosh 65649bd79ccSDebojyoti Ghosh /* --------------------------------------------------------------------------------------*/ 65749bd79ccSDebojyoti Ghosh 65849bd79ccSDebojyoti Ghosh /* zz = yy + Axx */ 659*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_SeqKAIJ(Mat A, Vec xx, Vec yy, Vec zz) 660d71ae5a4SJacob Faibussowitsch { 66149bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 66249bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; 66349bd79ccSDebojyoti Ghosh const PetscScalar *s = b->S, *t = b->T; 66449bd79ccSDebojyoti Ghosh const PetscScalar *x, *v, *bx; 66549bd79ccSDebojyoti Ghosh PetscScalar *y, *sums; 66649bd79ccSDebojyoti Ghosh const PetscInt m = b->AIJ->rmap->n, *idx, *ii; 66749bd79ccSDebojyoti Ghosh PetscInt n, i, jrow, j, l, p = b->p, q = b->q, k; 66849bd79ccSDebojyoti Ghosh 66949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 67049bd79ccSDebojyoti Ghosh if (!yy) { 6719566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 67249bd79ccSDebojyoti Ghosh } else { 6739566063dSJacob Faibussowitsch PetscCall(VecCopy(yy, zz)); 67449bd79ccSDebojyoti Ghosh } 6753ba16761SJacob Faibussowitsch if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(PETSC_SUCCESS); 67649bd79ccSDebojyoti Ghosh 6779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x)); 6789566063dSJacob Faibussowitsch PetscCall(VecGetArray(zz, &y)); 67949bd79ccSDebojyoti Ghosh idx = a->j; 68049bd79ccSDebojyoti Ghosh v = a->a; 68149bd79ccSDebojyoti Ghosh ii = a->i; 68249bd79ccSDebojyoti Ghosh 68349bd79ccSDebojyoti Ghosh if (b->isTI) { 68449bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 68549bd79ccSDebojyoti Ghosh jrow = ii[i]; 68649bd79ccSDebojyoti Ghosh n = ii[i + 1] - jrow; 68749bd79ccSDebojyoti Ghosh sums = y + p * i; 68849bd79ccSDebojyoti Ghosh for (j = 0; j < n; j++) { 689ad540459SPierre Jolivet for (k = 0; k < p; k++) sums[k] += v[jrow + j] * x[q * idx[jrow + j] + k]; 69049bd79ccSDebojyoti Ghosh } 69149bd79ccSDebojyoti Ghosh } 6929566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(3.0 * (a->nz) * p)); 69349bd79ccSDebojyoti Ghosh } else if (t) { 69449bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 69549bd79ccSDebojyoti Ghosh jrow = ii[i]; 69649bd79ccSDebojyoti Ghosh n = ii[i + 1] - jrow; 69749bd79ccSDebojyoti Ghosh sums = y + p * i; 69849bd79ccSDebojyoti Ghosh for (j = 0; j < n; j++) { 69949bd79ccSDebojyoti Ghosh for (k = 0; k < p; k++) { 700ad540459SPierre Jolivet for (l = 0; l < q; l++) sums[k] += v[jrow + j] * t[k + l * p] * x[q * idx[jrow + j] + l]; 70149bd79ccSDebojyoti Ghosh } 70249bd79ccSDebojyoti Ghosh } 70349bd79ccSDebojyoti Ghosh } 704234d9204SRichard Tran Mills /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do), 705234d9204SRichard Tran Mills * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T), 706234d9204SRichard Tran Mills * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter 707234d9204SRichard Tran Mills * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */ 7089566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((2.0 * p * q - p) * m + 2.0 * p * a->nz)); 70949bd79ccSDebojyoti Ghosh } 71049bd79ccSDebojyoti Ghosh if (s) { 71149bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 71249bd79ccSDebojyoti Ghosh sums = y + p * i; 71349bd79ccSDebojyoti Ghosh bx = x + q * i; 71449bd79ccSDebojyoti Ghosh if (i < b->AIJ->cmap->n) { 71549bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 716ad540459SPierre Jolivet for (k = 0; k < p; k++) sums[k] += s[k + j * p] * bx[j]; 71749bd79ccSDebojyoti Ghosh } 71849bd79ccSDebojyoti Ghosh } 71949bd79ccSDebojyoti Ghosh } 7209566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * m * p * q)); 72149bd79ccSDebojyoti Ghosh } 72249bd79ccSDebojyoti Ghosh 7239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x)); 7249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(zz, &y)); 7253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72649bd79ccSDebojyoti Ghosh } 72749bd79ccSDebojyoti Ghosh 728*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_SeqKAIJ(Mat A, Vec xx, Vec yy) 729d71ae5a4SJacob Faibussowitsch { 73049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 731f3fa974cSJacob Faibussowitsch PetscCall(MatMultAdd_SeqKAIJ(A, xx, NULL, yy)); 7323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73349bd79ccSDebojyoti Ghosh } 73449bd79ccSDebojyoti Ghosh 73549bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h> 73649bd79ccSDebojyoti Ghosh 737*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A, const PetscScalar **values) 738d71ae5a4SJacob Faibussowitsch { 73949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 74049bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; 74149bd79ccSDebojyoti Ghosh const PetscScalar *S = b->S; 74249bd79ccSDebojyoti Ghosh const PetscScalar *T = b->T; 74349bd79ccSDebojyoti Ghosh const PetscScalar *v = a->a; 74449bd79ccSDebojyoti Ghosh const PetscInt p = b->p, q = b->q, m = b->AIJ->rmap->n, *idx = a->j, *ii = a->i; 74549bd79ccSDebojyoti Ghosh PetscInt i, j, *v_pivots, dof, dof2; 74649bd79ccSDebojyoti Ghosh PetscScalar *diag, aval, *v_work; 74749bd79ccSDebojyoti Ghosh 74849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 74908401ef6SPierre Jolivet PetscCheck(p == q, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Block size must be square to calculate inverse."); 750aed4548fSBarry Smith PetscCheck(S || T || b->isTI, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATKAIJ: Cannot invert a zero matrix."); 75149bd79ccSDebojyoti Ghosh 75249bd79ccSDebojyoti Ghosh dof = p; 75349bd79ccSDebojyoti Ghosh dof2 = dof * dof; 75449bd79ccSDebojyoti Ghosh 75549bd79ccSDebojyoti Ghosh if (b->ibdiagvalid) { 75649bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 75849bd79ccSDebojyoti Ghosh } 7594dfa11a4SJacob Faibussowitsch if (!b->ibdiag) { PetscCall(PetscMalloc1(dof2 * m, &b->ibdiag)); } 76049bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 76149bd79ccSDebojyoti Ghosh diag = b->ibdiag; 76249bd79ccSDebojyoti Ghosh 7639566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dof, &v_work, dof, &v_pivots)); 76449bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 76549bd79ccSDebojyoti Ghosh if (S) { 7669566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(diag, S, dof2 * sizeof(PetscScalar))); 76749bd79ccSDebojyoti Ghosh } else { 7689566063dSJacob Faibussowitsch PetscCall(PetscMemzero(diag, dof2 * sizeof(PetscScalar))); 76949bd79ccSDebojyoti Ghosh } 77049bd79ccSDebojyoti Ghosh if (b->isTI) { 77149bd79ccSDebojyoti Ghosh aval = 0; 7729371c9d4SSatish Balay for (j = ii[i]; j < ii[i + 1]; j++) 7739371c9d4SSatish Balay if (idx[j] == i) aval = v[j]; 77449bd79ccSDebojyoti Ghosh for (j = 0; j < dof; j++) diag[j + dof * j] += aval; 77549bd79ccSDebojyoti Ghosh } else if (T) { 77649bd79ccSDebojyoti Ghosh aval = 0; 7779371c9d4SSatish Balay for (j = ii[i]; j < ii[i + 1]; j++) 7789371c9d4SSatish Balay if (idx[j] == i) aval = v[j]; 77949bd79ccSDebojyoti Ghosh for (j = 0; j < dof2; j++) diag[j] += aval * T[j]; 78049bd79ccSDebojyoti Ghosh } 7819566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A(dof, diag, v_pivots, v_work, PETSC_FALSE, NULL)); 78249bd79ccSDebojyoti Ghosh diag += dof2; 78349bd79ccSDebojyoti Ghosh } 7849566063dSJacob Faibussowitsch PetscCall(PetscFree2(v_work, v_pivots)); 78549bd79ccSDebojyoti Ghosh 78649bd79ccSDebojyoti Ghosh b->ibdiagvalid = PETSC_TRUE; 7873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78849bd79ccSDebojyoti Ghosh } 78949bd79ccSDebojyoti Ghosh 790d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A, Mat *B) 791d71ae5a4SJacob Faibussowitsch { 79249bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ *)A->data; 79349bd79ccSDebojyoti Ghosh 79449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 79549bd79ccSDebojyoti Ghosh *B = kaij->AIJ; 7963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79749bd79ccSDebojyoti Ghosh } 79849bd79ccSDebojyoti Ghosh 799d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatConvert_KAIJ_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) 800d71ae5a4SJacob Faibussowitsch { 801fac53fa1SPierre Jolivet Mat_SeqKAIJ *a = (Mat_SeqKAIJ *)A->data; 802fac53fa1SPierre Jolivet Mat AIJ, OAIJ, B; 803fac53fa1SPierre Jolivet PetscInt *d_nnz, *o_nnz = NULL, nz, i, j, m, d; 804fac53fa1SPierre Jolivet const PetscInt p = a->p, q = a->q; 805fac53fa1SPierre Jolivet PetscBool ismpikaij, missing; 806fac53fa1SPierre Jolivet 807fac53fa1SPierre Jolivet PetscFunctionBegin; 808fac53fa1SPierre Jolivet if (reuse != MAT_REUSE_MATRIX) { 8099566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIKAIJ, &ismpikaij)); 810fac53fa1SPierre Jolivet if (ismpikaij) { 811fac53fa1SPierre Jolivet Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 812fac53fa1SPierre Jolivet AIJ = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ; 813fac53fa1SPierre Jolivet OAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ; 814fac53fa1SPierre Jolivet } else { 815fac53fa1SPierre Jolivet AIJ = a->AIJ; 816fac53fa1SPierre Jolivet OAIJ = NULL; 817fac53fa1SPierre Jolivet } 8189566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 8199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 8209566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 8219566063dSJacob Faibussowitsch PetscCall(MatGetSize(AIJ, &m, NULL)); 8229566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(AIJ, &missing, &d)); /* assumption that all successive rows will have a missing diagonal */ 823fac53fa1SPierre Jolivet if (!missing || !a->S) d = m; 8249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &d_nnz)); 825fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8269566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(AIJ, i, &nz, NULL, NULL)); 827fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) d_nnz[i * p + j] = nz * q + (i >= d) * q; 8289566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(AIJ, i, &nz, NULL, NULL)); 829fac53fa1SPierre Jolivet } 830fac53fa1SPierre Jolivet if (OAIJ) { 8319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * p, &o_nnz)); 832fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8339566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL)); 834fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) o_nnz[i * p + j] = nz * q; 8359566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(OAIJ, i, &nz, NULL, NULL)); 836fac53fa1SPierre Jolivet } 8379566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B, 0, d_nnz, 0, o_nnz)); 838fac53fa1SPierre Jolivet } else { 8399566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(B, 0, d_nnz)); 840fac53fa1SPierre Jolivet } 8419566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz)); 8429566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nnz)); 843fac53fa1SPierre Jolivet } else B = *newmat; 8449566063dSJacob Faibussowitsch PetscCall(MatConvert_Basic(A, newtype, MAT_REUSE_MATRIX, &B)); 845fac53fa1SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 8469566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 847fac53fa1SPierre Jolivet } else *newmat = B; 8483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 849fac53fa1SPierre Jolivet } 850fac53fa1SPierre Jolivet 851*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatSOR_SeqKAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 852d71ae5a4SJacob Faibussowitsch { 85349bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *kaij = (Mat_SeqKAIJ *)A->data; 85449bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ *)kaij->AIJ->data; 85549bd79ccSDebojyoti Ghosh const PetscScalar *aa = a->a, *T = kaij->T, *v; 85649bd79ccSDebojyoti Ghosh const PetscInt m = kaij->AIJ->rmap->n, *ai = a->i, *aj = a->j, p = kaij->p, q = kaij->q, *diag, *vi; 85749bd79ccSDebojyoti Ghosh const PetscScalar *b, *xb, *idiag; 85849bd79ccSDebojyoti Ghosh PetscScalar *x, *work, *workt, *w, *y, *arr, *t, *arrt; 85949bd79ccSDebojyoti Ghosh PetscInt i, j, k, i2, bs, bs2, nz; 86049bd79ccSDebojyoti Ghosh 86149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 86249bd79ccSDebojyoti Ghosh its = its * lits; 863aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat"); 86408401ef6SPierre 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); 86528b400f6SJacob Faibussowitsch PetscCheck(!fshift, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for diagonal shift"); 86608401ef6SPierre 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"); 86708401ef6SPierre Jolivet PetscCheck(p == q, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatSOR for KAIJ: No support for non-square dense blocks"); 868f7d195e4SLawrence Mitchell bs = p; 869f7d195e4SLawrence Mitchell bs2 = bs * bs; 87049bd79ccSDebojyoti Ghosh 8713ba16761SJacob Faibussowitsch if (!m) PetscFunctionReturn(PETSC_SUCCESS); 87249bd79ccSDebojyoti Ghosh 8739566063dSJacob Faibussowitsch if (!kaij->ibdiagvalid) PetscCall(MatInvertBlockDiagonal_SeqKAIJ(A, NULL)); 87449bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 87549bd79ccSDebojyoti Ghosh diag = a->diag; 87649bd79ccSDebojyoti Ghosh 87749bd79ccSDebojyoti Ghosh if (!kaij->sor.setup) { 8789566063dSJacob 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)); 87949bd79ccSDebojyoti Ghosh kaij->sor.setup = PETSC_TRUE; 88049bd79ccSDebojyoti Ghosh } 88149bd79ccSDebojyoti Ghosh y = kaij->sor.y; 88249bd79ccSDebojyoti Ghosh w = kaij->sor.w; 88349bd79ccSDebojyoti Ghosh work = kaij->sor.work; 88449bd79ccSDebojyoti Ghosh t = kaij->sor.t; 88549bd79ccSDebojyoti Ghosh arr = kaij->sor.arr; 88649bd79ccSDebojyoti Ghosh 887d0609cedSBarry Smith PetscCall(VecGetArray(xx, &x)); 8889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 88949bd79ccSDebojyoti Ghosh 89049bd79ccSDebojyoti Ghosh if (flag & SOR_ZERO_INITIAL_GUESS) { 89149bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 89249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, b, idiag, x); /* x[0:bs] <- D^{-1} b[0:bs] */ 8939566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t, b, bs * sizeof(PetscScalar))); 89449bd79ccSDebojyoti Ghosh i2 = bs; 89549bd79ccSDebojyoti Ghosh idiag += bs2; 89649bd79ccSDebojyoti Ghosh for (i = 1; i < m; i++) { 89749bd79ccSDebojyoti Ghosh v = aa + ai[i]; 89849bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 89949bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 90049bd79ccSDebojyoti Ghosh 90149bd79ccSDebojyoti Ghosh if (T) { /* b - T (Arow * x) */ 9029566063dSJacob Faibussowitsch PetscCall(PetscMemzero(w, bs * sizeof(PetscScalar))); 90349bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 90449bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k]; 90549bd79ccSDebojyoti Ghosh } 90649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs, w, T, &t[i2]); 90749bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] += b[i2 + k]; 90849bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9099566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar))); 91049bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 91149bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) t[i2 + k] -= v[j] * x[vi[j] * bs + k]; 91249bd79ccSDebojyoti Ghosh } 91349bd79ccSDebojyoti Ghosh } else { 9149566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, b + i2, bs * sizeof(PetscScalar))); 91549bd79ccSDebojyoti Ghosh } 91649bd79ccSDebojyoti Ghosh 91749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, t + i2, idiag, y); 91849bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = omega * y[j]; 91949bd79ccSDebojyoti Ghosh 92049bd79ccSDebojyoti Ghosh idiag += bs2; 92149bd79ccSDebojyoti Ghosh i2 += bs; 92249bd79ccSDebojyoti Ghosh } 92349bd79ccSDebojyoti Ghosh /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 9249566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * a->nz)); 92549bd79ccSDebojyoti Ghosh xb = t; 92649bd79ccSDebojyoti Ghosh } else xb = b; 92749bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 92849bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1); 92949bd79ccSDebojyoti Ghosh i2 = bs * (m - 1); 9309566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 93149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, x + i2); 93249bd79ccSDebojyoti Ghosh i2 -= bs; 93349bd79ccSDebojyoti Ghosh idiag -= bs2; 93449bd79ccSDebojyoti Ghosh for (i = m - 2; i >= 0; i--) { 93549bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 93649bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 93749bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 93849bd79ccSDebojyoti Ghosh 93949bd79ccSDebojyoti Ghosh if (T) { /* FIXME: This branch untested */ 9409566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 94149bd79ccSDebojyoti Ghosh /* copy all rows of x that are needed into contiguous space */ 94249bd79ccSDebojyoti Ghosh workt = work; 94349bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9449566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 94549bd79ccSDebojyoti Ghosh workt += bs; 94649bd79ccSDebojyoti Ghosh } 94749bd79ccSDebojyoti Ghosh arrt = arr; 94849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9499566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 95049bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 95149bd79ccSDebojyoti Ghosh arrt += bs2; 95249bd79ccSDebojyoti Ghosh } 95349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 95449bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9559566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, t + i2, bs * sizeof(PetscScalar))); 95649bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 95749bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) w[k] -= v[j] * x[vi[j] * bs + k]; 95849bd79ccSDebojyoti Ghosh } 95949bd79ccSDebojyoti Ghosh } 96049bd79ccSDebojyoti Ghosh 96149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); /* RHS incorrect for omega != 1.0 */ 96249bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) x[i2 + j] = (1.0 - omega) * x[i2 + j] + omega * y[j]; 96349bd79ccSDebojyoti Ghosh 96449bd79ccSDebojyoti Ghosh idiag -= bs2; 96549bd79ccSDebojyoti Ghosh i2 -= bs; 96649bd79ccSDebojyoti Ghosh } 9679566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz))); 96849bd79ccSDebojyoti Ghosh } 96949bd79ccSDebojyoti Ghosh its--; 97049bd79ccSDebojyoti Ghosh } 97149bd79ccSDebojyoti Ghosh while (its--) { /* FIXME: This branch not updated */ 97249bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 97349bd79ccSDebojyoti Ghosh i2 = 0; 97449bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 97549bd79ccSDebojyoti Ghosh for (i = 0; i < m; i++) { 9769566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar))); 97749bd79ccSDebojyoti Ghosh 97849bd79ccSDebojyoti Ghosh v = aa + ai[i]; 97949bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 98049bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 98149bd79ccSDebojyoti Ghosh workt = work; 98249bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9839566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 98449bd79ccSDebojyoti Ghosh workt += bs; 98549bd79ccSDebojyoti Ghosh } 98649bd79ccSDebojyoti Ghosh arrt = arr; 98749bd79ccSDebojyoti Ghosh if (T) { 98849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9899566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 99049bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 99149bd79ccSDebojyoti Ghosh arrt += bs2; 99249bd79ccSDebojyoti Ghosh } 99349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 99449bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 99549bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 9969566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 99749bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 99849bd79ccSDebojyoti Ghosh arrt += bs2; 99949bd79ccSDebojyoti Ghosh } 100049bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 100149bd79ccSDebojyoti Ghosh } 10029566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(t + i2, w, bs * sizeof(PetscScalar))); 100349bd79ccSDebojyoti Ghosh 100449bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 100549bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 100649bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 100749bd79ccSDebojyoti Ghosh workt = work; 100849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10099566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 101049bd79ccSDebojyoti Ghosh workt += bs; 101149bd79ccSDebojyoti Ghosh } 101249bd79ccSDebojyoti Ghosh arrt = arr; 101349bd79ccSDebojyoti Ghosh if (T) { 101449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10159566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 101649bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 101749bd79ccSDebojyoti Ghosh arrt += bs2; 101849bd79ccSDebojyoti Ghosh } 101949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 102049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 102149bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10229566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 102349bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 102449bd79ccSDebojyoti Ghosh arrt += bs2; 102549bd79ccSDebojyoti Ghosh } 102649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 102749bd79ccSDebojyoti Ghosh } 102849bd79ccSDebojyoti Ghosh 102949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 103049bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 103149bd79ccSDebojyoti Ghosh 103249bd79ccSDebojyoti Ghosh idiag += bs2; 103349bd79ccSDebojyoti Ghosh i2 += bs; 103449bd79ccSDebojyoti Ghosh } 103549bd79ccSDebojyoti Ghosh xb = t; 10369371c9d4SSatish Balay } else xb = b; 103749bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 103849bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag + bs2 * (m - 1); 103949bd79ccSDebojyoti Ghosh i2 = bs * (m - 1); 104049bd79ccSDebojyoti Ghosh if (xb == b) { 104149bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) { 10429566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, b + i2, bs * sizeof(PetscScalar))); 104349bd79ccSDebojyoti Ghosh 104449bd79ccSDebojyoti Ghosh v = aa + ai[i]; 104549bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 104649bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 104749bd79ccSDebojyoti Ghosh workt = work; 104849bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10499566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 105049bd79ccSDebojyoti Ghosh workt += bs; 105149bd79ccSDebojyoti Ghosh } 105249bd79ccSDebojyoti Ghosh arrt = arr; 105349bd79ccSDebojyoti Ghosh if (T) { 105449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10559566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 105649bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 105749bd79ccSDebojyoti Ghosh arrt += bs2; 105849bd79ccSDebojyoti Ghosh } 105949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 106049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 106149bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10629566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 106349bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 106449bd79ccSDebojyoti Ghosh arrt += bs2; 106549bd79ccSDebojyoti Ghosh } 106649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 106749bd79ccSDebojyoti Ghosh } 106849bd79ccSDebojyoti Ghosh 106949bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 107049bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 107149bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 107249bd79ccSDebojyoti Ghosh workt = work; 107349bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10749566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 107549bd79ccSDebojyoti Ghosh workt += bs; 107649bd79ccSDebojyoti Ghosh } 107749bd79ccSDebojyoti Ghosh arrt = arr; 107849bd79ccSDebojyoti Ghosh if (T) { 107949bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10809566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 108149bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 108249bd79ccSDebojyoti Ghosh arrt += bs2; 108349bd79ccSDebojyoti Ghosh } 108449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 108549bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 108649bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 10879566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 108849bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 108949bd79ccSDebojyoti Ghosh arrt += bs2; 109049bd79ccSDebojyoti Ghosh } 109149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 109249bd79ccSDebojyoti Ghosh } 109349bd79ccSDebojyoti Ghosh 109449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 109549bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 109649bd79ccSDebojyoti Ghosh } 109749bd79ccSDebojyoti Ghosh } else { 109849bd79ccSDebojyoti Ghosh for (i = m - 1; i >= 0; i--) { 10999566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(w, xb + i2, bs * sizeof(PetscScalar))); 110049bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 110149bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 110249bd79ccSDebojyoti Ghosh nz = ai[i + 1] - diag[i] - 1; 110349bd79ccSDebojyoti Ghosh workt = work; 110449bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 11059566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(workt, x + bs * (*vi++), bs * sizeof(PetscScalar))); 110649bd79ccSDebojyoti Ghosh workt += bs; 110749bd79ccSDebojyoti Ghosh } 110849bd79ccSDebojyoti Ghosh arrt = arr; 110949bd79ccSDebojyoti Ghosh if (T) { 111049bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 11119566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(arrt, T, bs2 * sizeof(PetscScalar))); 111249bd79ccSDebojyoti Ghosh for (k = 0; k < bs2; k++) arrt[k] *= v[j]; 111349bd79ccSDebojyoti Ghosh arrt += bs2; 111449bd79ccSDebojyoti Ghosh } 111549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 111649bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 111749bd79ccSDebojyoti Ghosh for (j = 0; j < nz; j++) { 11189566063dSJacob Faibussowitsch PetscCall(PetscMemzero(arrt, bs2 * sizeof(PetscScalar))); 111949bd79ccSDebojyoti Ghosh for (k = 0; k < bs; k++) arrt[k + bs * k] = v[j]; 112049bd79ccSDebojyoti Ghosh arrt += bs2; 112149bd79ccSDebojyoti Ghosh } 112249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs, bs * nz, w, arr, work); 112349bd79ccSDebojyoti Ghosh } 112449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs, bs, w, idiag, y); 112549bd79ccSDebojyoti Ghosh for (j = 0; j < bs; j++) *(x + i2 + j) = (1.0 - omega) * *(x + i2 + j) + omega * *(y + j); 112649bd79ccSDebojyoti Ghosh } 112749bd79ccSDebojyoti Ghosh } 11289566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.0 * bs2 * (a->nz))); 112949bd79ccSDebojyoti Ghosh } 113049bd79ccSDebojyoti Ghosh } 113149bd79ccSDebojyoti Ghosh 1132d0609cedSBarry Smith PetscCall(VecRestoreArray(xx, &x)); 11339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 11343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113549bd79ccSDebojyoti Ghosh } 113649bd79ccSDebojyoti Ghosh 113749bd79ccSDebojyoti Ghosh /*===================================================================================*/ 113849bd79ccSDebojyoti Ghosh 1139*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMultAdd_MPIKAIJ(Mat A, Vec xx, Vec yy, Vec zz) 1140d71ae5a4SJacob Faibussowitsch { 114149bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 114249bd79ccSDebojyoti Ghosh 114349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 114449bd79ccSDebojyoti Ghosh if (!yy) { 11459566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 114649bd79ccSDebojyoti Ghosh } else { 11479566063dSJacob Faibussowitsch PetscCall(VecCopy(yy, zz)); 114849bd79ccSDebojyoti Ghosh } 11499566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */ 115049bd79ccSDebojyoti Ghosh /* start the scatter */ 11519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); 11529566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->multadd)(b->AIJ, xx, zz, zz)); 11539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); 11549566063dSJacob Faibussowitsch PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, zz, zz)); 11553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115649bd79ccSDebojyoti Ghosh } 115749bd79ccSDebojyoti Ghosh 1158*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatMult_MPIKAIJ(Mat A, Vec xx, Vec yy) 1159d71ae5a4SJacob Faibussowitsch { 116049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1161f3fa974cSJacob Faibussowitsch PetscCall(MatMultAdd_MPIKAIJ(A, xx, NULL, yy)); 11623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116349bd79ccSDebojyoti Ghosh } 116449bd79ccSDebojyoti Ghosh 1165*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A, const PetscScalar **values) 1166d71ae5a4SJacob Faibussowitsch { 116749bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 116849bd79ccSDebojyoti Ghosh 116949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 11709566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ is up to date. */ 11719566063dSJacob Faibussowitsch PetscCall((*b->AIJ->ops->invertblockdiagonal)(b->AIJ, values)); 11723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117349bd79ccSDebojyoti Ghosh } 117449bd79ccSDebojyoti Ghosh 117549bd79ccSDebojyoti Ghosh /* ----------------------------------------------------------------*/ 117649bd79ccSDebojyoti Ghosh 1177*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values) 1178d71ae5a4SJacob Faibussowitsch { 117949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ *)A->data; 11803ba16761SJacob Faibussowitsch PetscBool diag = PETSC_FALSE; 118149bd79ccSDebojyoti Ghosh PetscInt nzaij, nz, *colsaij, *idx, i, j, p = b->p, q = b->q, r = row / p, s = row % p, c; 118249bd79ccSDebojyoti Ghosh PetscScalar *vaij, *v, *S = b->S, *T = b->T; 118349bd79ccSDebojyoti Ghosh 118449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 118528b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active"); 118649bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 1187aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row); 118849bd79ccSDebojyoti Ghosh 118949bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 119049bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 119149bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 119249bd79ccSDebojyoti Ghosh if (values) *values = NULL; 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119449bd79ccSDebojyoti Ghosh } 119549bd79ccSDebojyoti Ghosh 119649bd79ccSDebojyoti Ghosh if (T || b->isTI) { 11979566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(b->AIJ, r, &nzaij, &colsaij, &vaij)); 119849bd79ccSDebojyoti Ghosh c = nzaij; 119949bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 120049bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 120149bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 120249bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 120349bd79ccSDebojyoti Ghosh c = i; 120449bd79ccSDebojyoti Ghosh } 120549bd79ccSDebojyoti Ghosh } 120649bd79ccSDebojyoti Ghosh } else nzaij = c = 0; 120749bd79ccSDebojyoti Ghosh 120849bd79ccSDebojyoti Ghosh /* calculate size of row */ 120949bd79ccSDebojyoti Ghosh nz = 0; 121049bd79ccSDebojyoti Ghosh if (S) nz += q; 121149bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (nzaij - 1) * q : nzaij * q); 121249bd79ccSDebojyoti Ghosh 121349bd79ccSDebojyoti Ghosh if (cols || values) { 12149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz, &idx, nz, &v)); 121538822f9dSRichard Tran Mills for (i = 0; i < q; i++) { 121638822f9dSRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 121738822f9dSRichard Tran Mills v[i] = 0.0; 121838822f9dSRichard Tran Mills } 121949bd79ccSDebojyoti Ghosh if (b->isTI) { 122049bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 122149bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 122249bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j; 122349bd79ccSDebojyoti Ghosh v[i * q + j] = (j == s ? vaij[i] : 0); 122449bd79ccSDebojyoti Ghosh } 122549bd79ccSDebojyoti Ghosh } 122649bd79ccSDebojyoti Ghosh } else if (T) { 122749bd79ccSDebojyoti Ghosh for (i = 0; i < nzaij; i++) { 122849bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 122949bd79ccSDebojyoti Ghosh idx[i * q + j] = colsaij[i] * q + j; 123049bd79ccSDebojyoti Ghosh v[i * q + j] = vaij[i] * T[s + j * p]; 123149bd79ccSDebojyoti Ghosh } 123249bd79ccSDebojyoti Ghosh } 123349bd79ccSDebojyoti Ghosh } 123449bd79ccSDebojyoti Ghosh if (S) { 123549bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 123649bd79ccSDebojyoti Ghosh idx[c * q + j] = r * q + j; 123749bd79ccSDebojyoti Ghosh v[c * q + j] += S[s + j * p]; 123849bd79ccSDebojyoti Ghosh } 123949bd79ccSDebojyoti Ghosh } 124049bd79ccSDebojyoti Ghosh } 124149bd79ccSDebojyoti Ghosh 124249bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 124349bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 124449bd79ccSDebojyoti Ghosh if (values) *values = v; 12453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124649bd79ccSDebojyoti Ghosh } 124749bd79ccSDebojyoti Ghosh 1248*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 1249d71ae5a4SJacob Faibussowitsch { 125049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1251cb4a9cd9SHong Zhang if (nz) *nz = 0; 12529566063dSJacob Faibussowitsch PetscCall(PetscFree2(*idx, *v)); 125349bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE; 12543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125549bd79ccSDebojyoti Ghosh } 125649bd79ccSDebojyoti Ghosh 1257*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *ncols, PetscInt **cols, PetscScalar **values) 1258d71ae5a4SJacob Faibussowitsch { 125949bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ *)A->data; 126049bd79ccSDebojyoti Ghosh Mat AIJ = b->A; 1261fc64b2cfSRichard Tran Mills PetscBool diag = PETSC_FALSE; 1262761d359dSRichard Tran Mills Mat MatAIJ, MatOAIJ; 126349bd79ccSDebojyoti Ghosh const PetscInt rstart = A->rmap->rstart, rend = A->rmap->rend, p = b->p, q = b->q, *garray; 1264389eba51SJed Brown PetscInt nz, *idx, ncolsaij = 0, ncolsoaij = 0, *colsaij, *colsoaij, r, s, c, i, j, lrow; 126549bd79ccSDebojyoti Ghosh PetscScalar *v, *vals, *ovals, *S = b->S, *T = b->T; 126649bd79ccSDebojyoti Ghosh 126749bd79ccSDebojyoti Ghosh PetscFunctionBegin; 12689566063dSJacob Faibussowitsch PetscCall(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */ 1269761d359dSRichard Tran Mills MatAIJ = ((Mat_SeqKAIJ *)b->AIJ->data)->AIJ; 1270761d359dSRichard Tran Mills MatOAIJ = ((Mat_SeqKAIJ *)b->OAIJ->data)->AIJ; 127128b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Already active"); 127249bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 1273aed4548fSBarry Smith PetscCheck(row >= rstart && row < rend, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only local rows"); 127449bd79ccSDebojyoti Ghosh lrow = row - rstart; 127549bd79ccSDebojyoti Ghosh 127649bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 127749bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 127849bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 127949bd79ccSDebojyoti Ghosh if (values) *values = NULL; 12803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 128149bd79ccSDebojyoti Ghosh } 128249bd79ccSDebojyoti Ghosh 128349bd79ccSDebojyoti Ghosh r = lrow / p; 128449bd79ccSDebojyoti Ghosh s = lrow % p; 128549bd79ccSDebojyoti Ghosh 128649bd79ccSDebojyoti Ghosh if (T || b->isTI) { 12879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJGetSeqAIJ(AIJ, NULL, NULL, &garray)); 12889566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(MatAIJ, lrow / p, &ncolsaij, &colsaij, &vals)); 12899566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(MatOAIJ, lrow / p, &ncolsoaij, &colsoaij, &ovals)); 129049bd79ccSDebojyoti Ghosh 129149bd79ccSDebojyoti Ghosh c = ncolsaij + ncolsoaij; 129249bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) { 129349bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 129449bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 129549bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 129649bd79ccSDebojyoti Ghosh c = i; 129749bd79ccSDebojyoti Ghosh } 129849bd79ccSDebojyoti Ghosh } 129949bd79ccSDebojyoti Ghosh } else c = 0; 130049bd79ccSDebojyoti Ghosh 130149bd79ccSDebojyoti Ghosh /* calculate size of row */ 130249bd79ccSDebojyoti Ghosh nz = 0; 130349bd79ccSDebojyoti Ghosh if (S) nz += q; 130449bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (ncolsaij + ncolsoaij - 1) * q : (ncolsaij + ncolsoaij) * q); 130549bd79ccSDebojyoti Ghosh 130649bd79ccSDebojyoti Ghosh if (cols || values) { 13079566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nz, &idx, nz, &v)); 1308a437a796SRichard Tran Mills for (i = 0; i < q; i++) { 1309a437a796SRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 1310a437a796SRichard Tran Mills v[i] = 0.0; 1311a437a796SRichard Tran Mills } 131249bd79ccSDebojyoti Ghosh if (b->isTI) { 131349bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) { 131449bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 131549bd79ccSDebojyoti Ghosh idx[i * q + j] = (colsaij[i] + rstart / p) * q + j; 131649bd79ccSDebojyoti Ghosh v[i * q + j] = (j == s ? vals[i] : 0.0); 131749bd79ccSDebojyoti Ghosh } 131849bd79ccSDebojyoti Ghosh } 131949bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsoaij; i++) { 132049bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 132149bd79ccSDebojyoti Ghosh idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j; 132249bd79ccSDebojyoti Ghosh v[(i + ncolsaij) * q + j] = (j == s ? ovals[i] : 0.0); 132349bd79ccSDebojyoti Ghosh } 132449bd79ccSDebojyoti Ghosh } 132549bd79ccSDebojyoti Ghosh } else if (T) { 132649bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsaij; i++) { 132749bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 132849bd79ccSDebojyoti Ghosh idx[i * q + j] = (colsaij[i] + rstart / p) * q + j; 132949bd79ccSDebojyoti Ghosh v[i * q + j] = vals[i] * T[s + j * p]; 133049bd79ccSDebojyoti Ghosh } 133149bd79ccSDebojyoti Ghosh } 133249bd79ccSDebojyoti Ghosh for (i = 0; i < ncolsoaij; i++) { 133349bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 133449bd79ccSDebojyoti Ghosh idx[(i + ncolsaij) * q + j] = garray[colsoaij[i]] * q + j; 133549bd79ccSDebojyoti Ghosh v[(i + ncolsaij) * q + j] = ovals[i] * T[s + j * p]; 133649bd79ccSDebojyoti Ghosh } 133749bd79ccSDebojyoti Ghosh } 133849bd79ccSDebojyoti Ghosh } 133949bd79ccSDebojyoti Ghosh if (S) { 134049bd79ccSDebojyoti Ghosh for (j = 0; j < q; j++) { 134149bd79ccSDebojyoti Ghosh idx[c * q + j] = (r + rstart / p) * q + j; 134249bd79ccSDebojyoti Ghosh v[c * q + j] += S[s + j * p]; 134349bd79ccSDebojyoti Ghosh } 134449bd79ccSDebojyoti Ghosh } 134549bd79ccSDebojyoti Ghosh } 134649bd79ccSDebojyoti Ghosh 134749bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 134849bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 134949bd79ccSDebojyoti Ghosh if (values) *values = v; 13503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135149bd79ccSDebojyoti Ghosh } 135249bd79ccSDebojyoti Ghosh 1353*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 1354d71ae5a4SJacob Faibussowitsch { 135549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 13569566063dSJacob Faibussowitsch PetscCall(PetscFree2(*idx, *v)); 135749bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ *)A->data)->getrowactive = PETSC_FALSE; 13583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135949bd79ccSDebojyoti Ghosh } 136049bd79ccSDebojyoti Ghosh 1361*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_KAIJ(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) 1362d71ae5a4SJacob Faibussowitsch { 136349bd79ccSDebojyoti Ghosh Mat A; 136449bd79ccSDebojyoti Ghosh 136549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 13669566063dSJacob Faibussowitsch PetscCall(MatConvert(mat, MATAIJ, MAT_INITIAL_MATRIX, &A)); 13679566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A, isrow, iscol, cll, newmat)); 13689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 13693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137049bd79ccSDebojyoti Ghosh } 137149bd79ccSDebojyoti Ghosh 137249bd79ccSDebojyoti Ghosh /* ---------------------------------------------------------------------------------- */ 137349bd79ccSDebojyoti Ghosh /*@C 137411a5261eSBarry Smith MatCreateKAIJ - Creates a matrix of type `MATKAIJ` to be used for matrices of the following form: 137549bd79ccSDebojyoti Ghosh 137649bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 137749bd79ccSDebojyoti Ghosh 137849bd79ccSDebojyoti Ghosh where 137949bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 138049bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 138111a5261eSBarry Smith A is an `MATAIJ` (n \times n) matrix 138249bd79ccSDebojyoti Ghosh I is the identity matrix 138349bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 138449bd79ccSDebojyoti Ghosh 138511a5261eSBarry Smith S and T are always stored independently on all processes as `PetscScalar` arrays in column-major format. 138649bd79ccSDebojyoti Ghosh 138749bd79ccSDebojyoti Ghosh Collective 138849bd79ccSDebojyoti Ghosh 138949bd79ccSDebojyoti Ghosh Input Parameters: 139011a5261eSBarry Smith + A - the `MATAIJ` matrix 139149bd79ccSDebojyoti Ghosh . p - number of rows in S and T 1392d3b90ce1SRichard Tran Mills . q - number of columns in S and T 139311a5261eSBarry Smith . S - the S matrix (can be NULL), stored as a `PetscScalar` array (column-major) 139411a5261eSBarry Smith - T - the T matrix (can be NULL), stored as a `PetscScalar` array (column-major) 139549bd79ccSDebojyoti Ghosh 139649bd79ccSDebojyoti Ghosh Output Parameter: 139711a5261eSBarry Smith . kaij - the new `MATKAIJ` matrix 139849bd79ccSDebojyoti Ghosh 1399d3b90ce1SRichard Tran Mills Notes: 140011a5261eSBarry 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. 140149bd79ccSDebojyoti Ghosh 140211a5261eSBarry Smith Changes to the entries of the `MATAIJ` matrix will immediately affect the `MATKAIJ` matrix. 140311a5261eSBarry Smith 140411a5261eSBarry Smith Developer Note: 140511a5261eSBarry Smith In the `MATMPIKAIJ` case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state 1406761d359dSRichard Tran Mills of the AIJ matrix 'A' that describes the blockwise action of the MATMPIKAIJ matrix and, if the object state has changed, lazily 1407761d359dSRichard Tran Mills rebuilding 'AIJ' and 'OAIJ' just before executing operations with the MATMPIKAIJ matrix. If new types of operations are added, 1408761d359dSRichard Tran Mills routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine). 1409761d359dSRichard Tran Mills 141049bd79ccSDebojyoti Ghosh Level: advanced 141149bd79ccSDebojyoti Ghosh 1412db781477SPatrick Sanan .seealso: `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MATKAIJ` 141349bd79ccSDebojyoti Ghosh @*/ 1414d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateKAIJ(Mat A, PetscInt p, PetscInt q, const PetscScalar S[], const PetscScalar T[], Mat *kaij) 1415d71ae5a4SJacob Faibussowitsch { 141649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 14179566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), kaij)); 14189566063dSJacob Faibussowitsch PetscCall(MatSetType(*kaij, MATKAIJ)); 14199566063dSJacob Faibussowitsch PetscCall(MatKAIJSetAIJ(*kaij, A)); 14209566063dSJacob Faibussowitsch PetscCall(MatKAIJSetS(*kaij, p, q, S)); 14219566063dSJacob Faibussowitsch PetscCall(MatKAIJSetT(*kaij, p, q, T)); 14229566063dSJacob Faibussowitsch PetscCall(MatSetUp(*kaij)); 14233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14240567c835SRichard Tran Mills } 142549bd79ccSDebojyoti Ghosh 14260567c835SRichard Tran Mills /*MC 14275881e567SRichard Tran Mills MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form 14285881e567SRichard Tran Mills [I \otimes S + A \otimes T], 14290567c835SRichard Tran Mills where 14305881e567SRichard Tran Mills S is a dense (p \times q) matrix, 14315881e567SRichard Tran Mills T is a dense (p \times q) matrix, 14325881e567SRichard Tran Mills A is an AIJ (n \times n) matrix, 14335881e567SRichard Tran Mills and I is the identity matrix. 14345881e567SRichard Tran Mills The resulting matrix is (np \times nq). 14350567c835SRichard Tran Mills 143611a5261eSBarry Smith S and T are always stored independently on all processes as `PetscScalar` arrays in column-major format. 14370567c835SRichard Tran Mills 143811a5261eSBarry Smith Note: 14395881e567SRichard 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, 14405881e567SRichard Tran Mills where x and b are column vectors containing the row-major representations of X and B. 14415881e567SRichard Tran Mills 14420567c835SRichard Tran Mills Level: advanced 14430567c835SRichard Tran Mills 1444db781477SPatrick Sanan .seealso: `MatKAIJSetAIJ()`, `MatKAIJSetS()`, `MatKAIJSetT()`, `MatKAIJGetAIJ()`, `MatKAIJGetS()`, `MatKAIJGetT()`, `MatCreateKAIJ()` 14450567c835SRichard Tran Mills M*/ 14460567c835SRichard Tran Mills 1447d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A) 1448d71ae5a4SJacob Faibussowitsch { 14490567c835SRichard Tran Mills Mat_MPIKAIJ *b; 14500567c835SRichard Tran Mills PetscMPIInt size; 14510567c835SRichard Tran Mills 14520567c835SRichard Tran Mills PetscFunctionBegin; 14534dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 14540567c835SRichard Tran Mills A->data = (void *)b; 14550567c835SRichard Tran Mills 14569566063dSJacob Faibussowitsch PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps))); 14570567c835SRichard Tran Mills 1458f4259b30SLisandro Dalcin b->w = NULL; 14599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 14600567c835SRichard Tran Mills if (size == 1) { 14619566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQKAIJ)); 14620567c835SRichard Tran Mills A->ops->destroy = MatDestroy_SeqKAIJ; 1463bb6fb833SRichard Tran Mills A->ops->mult = MatMult_SeqKAIJ; 1464bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_SeqKAIJ; 1465bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ; 14660567c835SRichard Tran Mills A->ops->getrow = MatGetRow_SeqKAIJ; 14670567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_SeqKAIJ; 14680567c835SRichard Tran Mills A->ops->sor = MatSOR_SeqKAIJ; 14699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqkaij_seqaij_C", MatConvert_KAIJ_AIJ)); 14700567c835SRichard Tran Mills } else { 14719566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATMPIKAIJ)); 14720567c835SRichard Tran Mills A->ops->destroy = MatDestroy_MPIKAIJ; 1473bb6fb833SRichard Tran Mills A->ops->mult = MatMult_MPIKAIJ; 1474bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_MPIKAIJ; 1475bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ; 14760567c835SRichard Tran Mills A->ops->getrow = MatGetRow_MPIKAIJ; 14770567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_MPIKAIJ; 14789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatGetDiagonalBlock_C", MatGetDiagonalBlock_MPIKAIJ)); 14799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_mpikaij_mpiaij_C", MatConvert_KAIJ_AIJ)); 14800567c835SRichard Tran Mills } 148106d61a37SPierre Jolivet A->ops->setup = MatSetUp_KAIJ; 148206d61a37SPierre Jolivet A->ops->view = MatView_KAIJ; 14830567c835SRichard Tran Mills A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ; 14843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148549bd79ccSDebojyoti Ghosh } 1486