149bd79ccSDebojyoti Ghosh 249bd79ccSDebojyoti Ghosh /* 349bd79ccSDebojyoti Ghosh Defines the basic matrix operations for the KAIJ matrix storage format. 449bd79ccSDebojyoti Ghosh This format is used to evaluate matrices of the form: 549bd79ccSDebojyoti Ghosh 649bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 749bd79ccSDebojyoti Ghosh 849bd79ccSDebojyoti Ghosh where 949bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 1049bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 1149bd79ccSDebojyoti Ghosh A is an AIJ (n \times n) matrix 1249bd79ccSDebojyoti Ghosh I is the identity matrix 1349bd79ccSDebojyoti Ghosh 1449bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 1549bd79ccSDebojyoti Ghosh 1649bd79ccSDebojyoti Ghosh We provide: 1749bd79ccSDebojyoti Ghosh MatMult() 1849bd79ccSDebojyoti Ghosh MatMultAdd() 1949bd79ccSDebojyoti Ghosh MatInvertBlockDiagonal() 2049bd79ccSDebojyoti Ghosh and 2149bd79ccSDebojyoti Ghosh MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*) 2249bd79ccSDebojyoti Ghosh 2349bd79ccSDebojyoti Ghosh This single directory handles both the sequential and parallel codes 2449bd79ccSDebojyoti Ghosh */ 2549bd79ccSDebojyoti Ghosh 2649bd79ccSDebojyoti Ghosh #include <../src/mat/impls/kaij/kaij.h> /*I "petscmat.h" I*/ 2749bd79ccSDebojyoti Ghosh #include <../src/mat/utils/freespace.h> 2849bd79ccSDebojyoti Ghosh #include <petsc/private/vecimpl.h> 2949bd79ccSDebojyoti Ghosh 3049bd79ccSDebojyoti Ghosh /*@C 3149bd79ccSDebojyoti Ghosh MatKAIJGetAIJ - Get the AIJ matrix describing the blockwise action of the KAIJ matrix 3249bd79ccSDebojyoti Ghosh 3349bd79ccSDebojyoti Ghosh Not Collective, but if the KAIJ matrix is parallel, the AIJ matrix is also parallel 3449bd79ccSDebojyoti Ghosh 3549bd79ccSDebojyoti Ghosh Input Parameter: 3649bd79ccSDebojyoti Ghosh . A - the KAIJ matrix 3749bd79ccSDebojyoti Ghosh 3849bd79ccSDebojyoti Ghosh Output Parameter: 3949bd79ccSDebojyoti Ghosh . B - the AIJ matrix 4049bd79ccSDebojyoti Ghosh 4149bd79ccSDebojyoti Ghosh Level: advanced 4249bd79ccSDebojyoti Ghosh 4349bd79ccSDebojyoti Ghosh Notes: The reference count on the AIJ matrix is not increased so you should not destroy it. 4449bd79ccSDebojyoti Ghosh 4549bd79ccSDebojyoti Ghosh .seealso: MatCreateKAIJ() 4649bd79ccSDebojyoti Ghosh @*/ 4749bd79ccSDebojyoti Ghosh PetscErrorCode MatKAIJGetAIJ(Mat A,Mat *B) 4849bd79ccSDebojyoti Ghosh { 4949bd79ccSDebojyoti Ghosh PetscBool ismpikaij,isseqkaij; 5049bd79ccSDebojyoti Ghosh 5149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQKAIJ,&isseqkaij)); 5449bd79ccSDebojyoti Ghosh if (ismpikaij) { 5549bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 5649bd79ccSDebojyoti Ghosh 5749bd79ccSDebojyoti Ghosh *B = b->A; 5849bd79ccSDebojyoti Ghosh } else if (isseqkaij) { 5949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 6049bd79ccSDebojyoti Ghosh 6149bd79ccSDebojyoti Ghosh *B = b->AIJ; 62b04351cbSRichard Tran Mills } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix passed in is not of type KAIJ"); 6349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 6449bd79ccSDebojyoti Ghosh } 6549bd79ccSDebojyoti Ghosh 6649bd79ccSDebojyoti Ghosh /*@C 6749bd79ccSDebojyoti Ghosh MatKAIJGetS - Get the S matrix describing the shift action of the KAIJ matrix 6849bd79ccSDebojyoti Ghosh 690567c835SRichard Tran Mills Not Collective; the entire S is stored and returned independently on all processes. 7049bd79ccSDebojyoti Ghosh 7149bd79ccSDebojyoti Ghosh Input Parameter: 7249bd79ccSDebojyoti Ghosh . A - the KAIJ matrix 7349bd79ccSDebojyoti Ghosh 74a5b5c723SRichard Tran Mills Output Parameters: 75a5b5c723SRichard Tran Mills + m - the number of rows in S 76a5b5c723SRichard Tran Mills . n - the number of columns in S 77a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 7849bd79ccSDebojyoti Ghosh 79a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 8031a97b9aSRichard Tran Mills 8149bd79ccSDebojyoti Ghosh Level: advanced 8249bd79ccSDebojyoti Ghosh 8331a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 8449bd79ccSDebojyoti Ghosh @*/ 85a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetS(Mat A,PetscInt *m,PetscInt *n,PetscScalar **S) 8649bd79ccSDebojyoti Ghosh { 8749bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 8849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 89a5b5c723SRichard Tran Mills if (m) *m = b->p; 90a5b5c723SRichard Tran Mills if (n) *n = b->q; 91a5b5c723SRichard Tran Mills if (S) *S = b->S; 92a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 93a5b5c723SRichard Tran Mills } 94a5b5c723SRichard Tran Mills 95a5b5c723SRichard Tran Mills /*@C 96a5b5c723SRichard Tran Mills MatKAIJGetSRead - Get a read-only pointer to the S matrix describing the shift action of the KAIJ matrix 97a5b5c723SRichard Tran Mills 98a5b5c723SRichard Tran Mills Not Collective; the entire S is stored and returned independently on all processes. 99a5b5c723SRichard Tran Mills 100a5b5c723SRichard Tran Mills Input Parameter: 101a5b5c723SRichard Tran Mills . A - the KAIJ matrix 102a5b5c723SRichard Tran Mills 103a5b5c723SRichard Tran Mills Output Parameters: 104a5b5c723SRichard Tran Mills + m - the number of rows in S 105a5b5c723SRichard Tran Mills . n - the number of columns in S 106a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 107a5b5c723SRichard Tran Mills 108a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 109a5b5c723SRichard Tran Mills 110a5b5c723SRichard Tran Mills Level: advanced 111a5b5c723SRichard Tran Mills 112a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 113a5b5c723SRichard Tran Mills @*/ 114a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetSRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **S) 115a5b5c723SRichard Tran Mills { 116a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 117a5b5c723SRichard Tran Mills PetscFunctionBegin; 118a5b5c723SRichard Tran Mills if (m) *m = b->p; 119a5b5c723SRichard Tran Mills if (n) *n = b->q; 120a5b5c723SRichard Tran Mills if (S) *S = b->S; 121a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 122a5b5c723SRichard Tran Mills } 123a5b5c723SRichard Tran Mills 124a5b5c723SRichard Tran Mills /*@C 125a5b5c723SRichard Tran Mills MatKAIJRestoreS - Restore array obtained with MatKAIJGetS() 126a5b5c723SRichard Tran Mills 127a5b5c723SRichard Tran Mills Not collective 128a5b5c723SRichard Tran Mills 129a5b5c723SRichard Tran Mills Input Parameter: 130a5b5c723SRichard Tran Mills . A - the KAIJ matrix 131a5b5c723SRichard Tran Mills 132a5b5c723SRichard Tran Mills Output Parameter: 133a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS() 134a5b5c723SRichard Tran Mills 135a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 136a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 137a5b5c723SRichard Tran Mills 138a5b5c723SRichard Tran Mills Level: advanced 139a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead() 140a5b5c723SRichard Tran Mills @*/ 141a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreS(Mat A,PetscScalar **S) 142a5b5c723SRichard Tran Mills { 143a5b5c723SRichard Tran Mills PetscFunctionBegin; 144a5b5c723SRichard Tran Mills if (S) *S = NULL; 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateIncrease((PetscObject)A)); 146a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 147a5b5c723SRichard Tran Mills } 148a5b5c723SRichard Tran Mills 149a5b5c723SRichard Tran Mills /*@C 150a5b5c723SRichard Tran Mills MatKAIJRestoreSRead - Restore array obtained with MatKAIJGetSRead() 151a5b5c723SRichard Tran Mills 152a5b5c723SRichard Tran Mills Not collective 153a5b5c723SRichard Tran Mills 154a5b5c723SRichard Tran Mills Input Parameter: 155a5b5c723SRichard Tran Mills . A - the KAIJ matrix 156a5b5c723SRichard Tran Mills 157a5b5c723SRichard Tran Mills Output Parameter: 158a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS() 159a5b5c723SRichard Tran Mills 160a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 161a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 162a5b5c723SRichard Tran Mills 163a5b5c723SRichard Tran Mills Level: advanced 164a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead() 165a5b5c723SRichard Tran Mills @*/ 166a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreSRead(Mat A,const PetscScalar **S) 167a5b5c723SRichard Tran Mills { 168a5b5c723SRichard Tran Mills PetscFunctionBegin; 169a5b5c723SRichard Tran Mills if (S) *S = NULL; 17049bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 17149bd79ccSDebojyoti Ghosh } 17249bd79ccSDebojyoti Ghosh 17349bd79ccSDebojyoti Ghosh /*@C 17431a97b9aSRichard Tran Mills MatKAIJGetT - Get the transformation matrix T associated with the KAIJ matrix 17549bd79ccSDebojyoti Ghosh 1760567c835SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 17749bd79ccSDebojyoti Ghosh 17849bd79ccSDebojyoti Ghosh Input Parameter: 17949bd79ccSDebojyoti Ghosh . A - the KAIJ matrix 18049bd79ccSDebojyoti Ghosh 181d8d19677SJose E. Roman Output Parameters: 182a5b5c723SRichard Tran Mills + m - the number of rows in T 183a5b5c723SRichard Tran Mills . n - the number of columns in T 184a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 18549bd79ccSDebojyoti Ghosh 186a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 18731a97b9aSRichard Tran Mills 18849bd79ccSDebojyoti Ghosh Level: advanced 18949bd79ccSDebojyoti Ghosh 19031a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 19149bd79ccSDebojyoti Ghosh @*/ 192a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetT(Mat A,PetscInt *m,PetscInt *n,PetscScalar **T) 19349bd79ccSDebojyoti Ghosh { 19449bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 19549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 196a5b5c723SRichard Tran Mills if (m) *m = b->p; 197a5b5c723SRichard Tran Mills if (n) *n = b->q; 198a5b5c723SRichard Tran Mills if (T) *T = b->T; 199a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 200a5b5c723SRichard Tran Mills } 201a5b5c723SRichard Tran Mills 202a5b5c723SRichard Tran Mills /*@C 203a5b5c723SRichard Tran Mills MatKAIJGetTRead - Get a read-only pointer to the transformation matrix T associated with the KAIJ matrix 204a5b5c723SRichard Tran Mills 205a5b5c723SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 206a5b5c723SRichard Tran Mills 207a5b5c723SRichard Tran Mills Input Parameter: 208a5b5c723SRichard Tran Mills . A - the KAIJ matrix 209a5b5c723SRichard Tran Mills 210d8d19677SJose E. Roman Output Parameters: 211a5b5c723SRichard Tran Mills + m - the number of rows in T 212a5b5c723SRichard Tran Mills . n - the number of columns in T 213a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 214a5b5c723SRichard Tran Mills 215a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 216a5b5c723SRichard Tran Mills 217a5b5c723SRichard Tran Mills Level: advanced 218a5b5c723SRichard Tran Mills 219a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 220a5b5c723SRichard Tran Mills @*/ 221a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetTRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **T) 222a5b5c723SRichard Tran Mills { 223a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 224a5b5c723SRichard Tran Mills PetscFunctionBegin; 225a5b5c723SRichard Tran Mills if (m) *m = b->p; 226a5b5c723SRichard Tran Mills if (n) *n = b->q; 227a5b5c723SRichard Tran Mills if (T) *T = b->T; 228a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 229a5b5c723SRichard Tran Mills } 230a5b5c723SRichard Tran Mills 231a5b5c723SRichard Tran Mills /*@C 232a5b5c723SRichard Tran Mills MatKAIJRestoreT - Restore array obtained with MatKAIJGetT() 233a5b5c723SRichard Tran Mills 234a5b5c723SRichard Tran Mills Not collective 235a5b5c723SRichard Tran Mills 236a5b5c723SRichard Tran Mills Input Parameter: 237a5b5c723SRichard Tran Mills . A - the KAIJ matrix 238a5b5c723SRichard Tran Mills 239a5b5c723SRichard Tran Mills Output Parameter: 240a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS() 241a5b5c723SRichard Tran Mills 242a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 243a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 244a5b5c723SRichard Tran Mills 245a5b5c723SRichard Tran Mills Level: advanced 246a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead() 247a5b5c723SRichard Tran Mills @*/ 248a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreT(Mat A,PetscScalar **T) 249a5b5c723SRichard Tran Mills { 250a5b5c723SRichard Tran Mills PetscFunctionBegin; 251a5b5c723SRichard Tran Mills if (T) *T = NULL; 2525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateIncrease((PetscObject)A)); 253a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 254a5b5c723SRichard Tran Mills } 255a5b5c723SRichard Tran Mills 256a5b5c723SRichard Tran Mills /*@C 257a5b5c723SRichard Tran Mills MatKAIJRestoreTRead - Restore array obtained with MatKAIJGetTRead() 258a5b5c723SRichard Tran Mills 259a5b5c723SRichard Tran Mills Not collective 260a5b5c723SRichard Tran Mills 261a5b5c723SRichard Tran Mills Input Parameter: 262a5b5c723SRichard Tran Mills . A - the KAIJ matrix 263a5b5c723SRichard Tran Mills 264a5b5c723SRichard Tran Mills Output Parameter: 265a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS() 266a5b5c723SRichard Tran Mills 267a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 268a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 269a5b5c723SRichard Tran Mills 270a5b5c723SRichard Tran Mills Level: advanced 271a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead() 272a5b5c723SRichard Tran Mills @*/ 273a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreTRead(Mat A,const PetscScalar **T) 274a5b5c723SRichard Tran Mills { 275a5b5c723SRichard Tran Mills PetscFunctionBegin; 276a5b5c723SRichard Tran Mills if (T) *T = NULL; 27749bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 27849bd79ccSDebojyoti Ghosh } 27949bd79ccSDebojyoti Ghosh 2800567c835SRichard Tran Mills /*@ 2810567c835SRichard Tran Mills MatKAIJSetAIJ - Set the AIJ matrix describing the blockwise action of the KAIJ matrix 2820567c835SRichard Tran Mills 2830567c835SRichard Tran Mills Logically Collective; if the AIJ matrix is parallel, the KAIJ matrix is also parallel 2840567c835SRichard Tran Mills 2850567c835SRichard Tran Mills Input Parameters: 2860567c835SRichard Tran Mills + A - the KAIJ matrix 2870567c835SRichard Tran Mills - B - the AIJ matrix 2880567c835SRichard Tran Mills 28915b9d025SRichard Tran Mills Notes: 29015b9d025SRichard Tran Mills This function increases the reference count on the AIJ matrix, so the user is free to destroy the matrix if it is not needed. 29115b9d025SRichard Tran Mills Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix. 29215b9d025SRichard Tran Mills 2930567c835SRichard Tran Mills Level: advanced 2940567c835SRichard Tran Mills 2950567c835SRichard Tran Mills .seealso: MatKAIJGetAIJ(), MatKAIJSetS(), MatKAIJSetT() 2960567c835SRichard Tran Mills @*/ 2970567c835SRichard Tran Mills PetscErrorCode MatKAIJSetAIJ(Mat A,Mat B) 2980567c835SRichard Tran Mills { 2990567c835SRichard Tran Mills PetscMPIInt size; 30006d61a37SPierre Jolivet PetscBool flg; 3010567c835SRichard Tran Mills 3020567c835SRichard Tran Mills PetscFunctionBegin; 3035f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 3040567c835SRichard Tran Mills if (size == 1) { 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&flg)); 306*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatKAIJSetAIJ() with MATSEQKAIJ does not support %s as the AIJ mat",((PetscObject)B)->type_name); 3070567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 3080567c835SRichard Tran Mills a->AIJ = B; 3090567c835SRichard Tran Mills } else { 3100567c835SRichard Tran Mills Mat_MPIKAIJ *a = (Mat_MPIKAIJ*)A->data; 3110567c835SRichard Tran Mills a->A = B; 3120567c835SRichard Tran Mills } 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)B)); 3140567c835SRichard Tran Mills PetscFunctionReturn(0); 3150567c835SRichard Tran Mills } 3160567c835SRichard Tran Mills 3170567c835SRichard Tran Mills /*@C 3180567c835SRichard Tran Mills MatKAIJSetS - Set the S matrix describing the shift action of the KAIJ matrix 3190567c835SRichard Tran Mills 3200567c835SRichard Tran Mills Logically Collective; the entire S is stored independently on all processes. 3210567c835SRichard Tran Mills 3220567c835SRichard Tran Mills Input Parameters: 3230567c835SRichard Tran Mills + A - the KAIJ matrix 3240567c835SRichard Tran Mills . p - the number of rows in S 3250567c835SRichard Tran Mills . q - the number of columns in S 3260567c835SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 3270567c835SRichard Tran Mills 3280567c835SRichard Tran Mills Notes: The dimensions p and q must match those of the transformation matrix T associated with the KAIJ matrix. 32988f48298SRichard Tran Mills The S matrix is copied, so the user can destroy this array. 3300567c835SRichard Tran Mills 3310567c835SRichard Tran Mills Level: Advanced 3320567c835SRichard Tran Mills 3330567c835SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJSetT(), MatKAIJSetAIJ() 3340567c835SRichard Tran Mills @*/ 3350567c835SRichard Tran Mills PetscErrorCode MatKAIJSetS(Mat A,PetscInt p,PetscInt q,const PetscScalar S[]) 3360567c835SRichard Tran Mills { 3370567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 3380567c835SRichard Tran Mills 3390567c835SRichard Tran Mills PetscFunctionBegin; 3405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(a->S)); 3410567c835SRichard Tran Mills if (S) { 3425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(p*q,&a->S)); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(a->S,S,p*q*sizeof(PetscScalar))); 3440567c835SRichard Tran Mills } else a->S = NULL; 3450567c835SRichard Tran Mills 3460567c835SRichard Tran Mills a->p = p; 3470567c835SRichard Tran Mills a->q = q; 3480567c835SRichard Tran Mills PetscFunctionReturn(0); 3490567c835SRichard Tran Mills } 3500567c835SRichard Tran Mills 3510567c835SRichard Tran Mills /*@C 352910cf402Sprj- MatKAIJGetScaledIdentity - Check if both S and T are scaled identities. 353910cf402Sprj- 354910cf402Sprj- Logically Collective. 355910cf402Sprj- 356910cf402Sprj- Input Parameter: 357910cf402Sprj- . A - the KAIJ matrix 358910cf402Sprj- 359910cf402Sprj- Output Parameter: 360910cf402Sprj- . identity - the Boolean value 361910cf402Sprj- 362910cf402Sprj- Level: Advanced 363910cf402Sprj- 364910cf402Sprj- .seealso: MatKAIJGetS(), MatKAIJGetT() 365910cf402Sprj- @*/ 366910cf402Sprj- PetscErrorCode MatKAIJGetScaledIdentity(Mat A,PetscBool* identity) 367910cf402Sprj- { 368910cf402Sprj- Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 369910cf402Sprj- PetscInt i,j; 370910cf402Sprj- 371910cf402Sprj- PetscFunctionBegin; 372910cf402Sprj- if (a->p != a->q) { 373910cf402Sprj- *identity = PETSC_FALSE; 374910cf402Sprj- PetscFunctionReturn(0); 375910cf402Sprj- } else *identity = PETSC_TRUE; 376910cf402Sprj- if (!a->isTI || a->S) { 377910cf402Sprj- for (i=0; i<a->p && *identity; i++) { 378910cf402Sprj- for (j=0; j<a->p && *identity; j++) { 379910cf402Sprj- if (i != j) { 380910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 381910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 382910cf402Sprj- } else { 383910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i*(a->p+1)]-a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 384910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i*(a->p+1)]-a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 385910cf402Sprj- } 386910cf402Sprj- } 387910cf402Sprj- } 388910cf402Sprj- } 389910cf402Sprj- PetscFunctionReturn(0); 390910cf402Sprj- } 391910cf402Sprj- 392910cf402Sprj- /*@C 3930567c835SRichard Tran Mills MatKAIJSetT - Set the transformation matrix T associated with the KAIJ matrix 3940567c835SRichard Tran Mills 3950567c835SRichard Tran Mills Logically Collective; the entire T is stored independently on all processes. 3960567c835SRichard Tran Mills 3970567c835SRichard Tran Mills Input Parameters: 3980567c835SRichard Tran Mills + A - the KAIJ matrix 3990567c835SRichard Tran Mills . p - the number of rows in S 4000567c835SRichard Tran Mills . q - the number of columns in S 4010567c835SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 4020567c835SRichard Tran Mills 4030567c835SRichard Tran Mills Notes: The dimensions p and q must match those of the shift matrix S associated with the KAIJ matrix. 40488f48298SRichard Tran Mills The T matrix is copied, so the user can destroy this array. 4050567c835SRichard Tran Mills 4060567c835SRichard Tran Mills Level: Advanced 4070567c835SRichard Tran Mills 4080567c835SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJSetS(), MatKAIJSetAIJ() 4090567c835SRichard Tran Mills @*/ 4100567c835SRichard Tran Mills PetscErrorCode MatKAIJSetT(Mat A,PetscInt p,PetscInt q,const PetscScalar T[]) 4110567c835SRichard Tran Mills { 4120567c835SRichard Tran Mills PetscInt i,j; 4130567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 4140567c835SRichard Tran Mills PetscBool isTI = PETSC_FALSE; 4150567c835SRichard Tran Mills 4160567c835SRichard Tran Mills PetscFunctionBegin; 4170567c835SRichard Tran Mills /* check if T is an identity matrix */ 4180567c835SRichard Tran Mills if (T && (p == q)) { 4190567c835SRichard Tran Mills isTI = PETSC_TRUE; 4200567c835SRichard Tran Mills for (i=0; i<p; i++) { 4210567c835SRichard Tran Mills for (j=0; j<q; j++) { 4220567c835SRichard Tran Mills if (i == j) { 4230567c835SRichard Tran Mills /* diagonal term must be 1 */ 4240567c835SRichard Tran Mills if (T[i+j*p] != 1.0) isTI = PETSC_FALSE; 4250567c835SRichard Tran Mills } else { 4260567c835SRichard Tran Mills /* off-diagonal term must be 0 */ 4270567c835SRichard Tran Mills if (T[i+j*p] != 0.0) isTI = PETSC_FALSE; 4280567c835SRichard Tran Mills } 4290567c835SRichard Tran Mills } 4300567c835SRichard Tran Mills } 4310567c835SRichard Tran Mills } 4320567c835SRichard Tran Mills a->isTI = isTI; 4330567c835SRichard Tran Mills 4345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(a->T)); 4350567c835SRichard Tran Mills if (T && (!isTI)) { 4365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(p*q,&a->T)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(a->T,T,p*q*sizeof(PetscScalar))); 43850d19d74SRichard Tran Mills } else a->T = NULL; 4390567c835SRichard Tran Mills 4400567c835SRichard Tran Mills a->p = p; 4410567c835SRichard Tran Mills a->q = q; 4420567c835SRichard Tran Mills PetscFunctionReturn(0); 4430567c835SRichard Tran Mills } 4440567c835SRichard Tran Mills 44549bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_SeqKAIJ(Mat A) 44649bd79ccSDebojyoti Ghosh { 44749bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 44849bd79ccSDebojyoti Ghosh 44949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 4505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&b->AIJ)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->S)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->T)); 4535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->ibdiag)); 4545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree5(b->sor.w,b->sor.y,b->sor.work,b->sor.t,b->sor.arr)); 4555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqkaij_seqaij_C",NULL)); 4565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(A->data)); 45749bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 45849bd79ccSDebojyoti Ghosh } 45949bd79ccSDebojyoti Ghosh 460e0e5a793SRichard Tran Mills PETSC_INTERN PetscErrorCode MatKAIJ_build_AIJ_OAIJ(Mat A) 461e0e5a793SRichard Tran Mills { 462e0e5a793SRichard Tran Mills Mat_MPIKAIJ *a; 463e0e5a793SRichard Tran Mills Mat_MPIAIJ *mpiaij; 464e0e5a793SRichard Tran Mills PetscScalar *T; 465e0e5a793SRichard Tran Mills PetscInt i,j; 466e0e5a793SRichard Tran Mills PetscObjectState state; 467e0e5a793SRichard Tran Mills 468e0e5a793SRichard Tran Mills PetscFunctionBegin; 469e0e5a793SRichard Tran Mills a = (Mat_MPIKAIJ*)A->data; 470e0e5a793SRichard Tran Mills mpiaij = (Mat_MPIAIJ*)a->A->data; 471e0e5a793SRichard Tran Mills 4725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectStateGet((PetscObject)a->A,&state)); 473e0e5a793SRichard Tran Mills if (state == a->state) { 474e0e5a793SRichard Tran Mills /* The existing AIJ and KAIJ members are up-to-date, so simply exit. */ 475e0e5a793SRichard Tran Mills PetscFunctionReturn(0); 476e0e5a793SRichard Tran Mills } else { 4775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&a->AIJ)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&a->OAIJ)); 479e0e5a793SRichard Tran Mills if (a->isTI) { 480e0e5a793SRichard Tran Mills /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL. 481e0e5a793SRichard Tran Mills * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will 482e0e5a793SRichard Tran Mills * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix 483e0e5a793SRichard Tran Mills * to pass in. */ 4845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(a->p*a->q,&T)); 485e0e5a793SRichard Tran Mills for (i=0; i<a->p; i++) { 486e0e5a793SRichard Tran Mills for (j=0; j<a->q; j++) { 487e0e5a793SRichard Tran Mills if (i==j) T[i+j*a->p] = 1.0; 488e0e5a793SRichard Tran Mills else T[i+j*a->p] = 0.0; 489e0e5a793SRichard Tran Mills } 490e0e5a793SRichard Tran Mills } 491e0e5a793SRichard Tran Mills } else T = a->T; 4925f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(mpiaij->A,a->p,a->q,a->S,T,&a->AIJ)); 4935f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(mpiaij->B,a->p,a->q,NULL,T,&a->OAIJ)); 494e0e5a793SRichard Tran Mills if (a->isTI) { 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(T)); 496e0e5a793SRichard Tran Mills } 497e0e5a793SRichard Tran Mills a->state = state; 498e0e5a793SRichard Tran Mills } 499e0e5a793SRichard Tran Mills 500e0e5a793SRichard Tran Mills PetscFunctionReturn(0); 501e0e5a793SRichard Tran Mills } 502e0e5a793SRichard Tran Mills 50349bd79ccSDebojyoti Ghosh PetscErrorCode MatSetUp_KAIJ(Mat A) 50449bd79ccSDebojyoti Ghosh { 5050567c835SRichard Tran Mills PetscInt n; 5060567c835SRichard Tran Mills PetscMPIInt size; 5070567c835SRichard Tran Mills Mat_SeqKAIJ *seqkaij = (Mat_SeqKAIJ*)A->data; 5080567c835SRichard Tran Mills 50949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5105f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 5110567c835SRichard Tran Mills if (size == 1) { 5125f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 5135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetBlockSize(A->rmap,seqkaij->p)); 5145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetBlockSize(A->cmap,seqkaij->q)); 5155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(A->rmap)); 5165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(A->cmap)); 5170567c835SRichard Tran Mills } else { 5180567c835SRichard Tran Mills Mat_MPIKAIJ *a; 5190567c835SRichard Tran Mills Mat_MPIAIJ *mpiaij; 5200567c835SRichard Tran Mills IS from,to; 5210567c835SRichard Tran Mills Vec gvec; 5220567c835SRichard Tran Mills 5230567c835SRichard Tran Mills a = (Mat_MPIKAIJ*)A->data; 524d3f912faSRichard Tran Mills mpiaij = (Mat_MPIAIJ*)a->A->data; 5255f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 5265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetBlockSize(A->rmap,seqkaij->p)); 5275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetBlockSize(A->cmap,seqkaij->q)); 5285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(A->rmap)); 5295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(A->cmap)); 5300567c835SRichard Tran Mills 5315f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJ_build_AIJ_OAIJ(A)); 5320567c835SRichard Tran Mills 5335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(mpiaij->lvec,&n)); 5345f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&a->w)); 5355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(a->w,n*a->q,n*a->q)); 5365f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(a->w,a->q)); 5375f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetType(a->w,VECSEQ)); 5380567c835SRichard Tran Mills 5390567c835SRichard Tran Mills /* create two temporary Index sets for build scatter gather */ 5405f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateBlock(PetscObjectComm((PetscObject)a->A),a->q,n,mpiaij->garray,PETSC_COPY_VALUES,&from)); 5415f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,n*a->q,0,1,&to)); 5420567c835SRichard Tran Mills 5430567c835SRichard Tran Mills /* create temporary global vector to generate scatter context */ 5445f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A),a->q,a->q*a->A->cmap->n,a->q*a->A->cmap->N,NULL,&gvec)); 5450567c835SRichard Tran Mills 5460567c835SRichard Tran Mills /* generate the scatter context */ 5475f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(gvec,from,a->w,to,&a->ctx)); 5480567c835SRichard Tran Mills 5495f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 5505f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 5515f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&gvec)); 5520567c835SRichard Tran Mills } 5530567c835SRichard Tran Mills 5540567c835SRichard Tran Mills A->assembled = PETSC_TRUE; 55549bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 55649bd79ccSDebojyoti Ghosh } 55749bd79ccSDebojyoti Ghosh 558e6985dafSRichard Tran Mills PetscErrorCode MatView_KAIJ(Mat A,PetscViewer viewer) 55949bd79ccSDebojyoti Ghosh { 560e6985dafSRichard Tran Mills PetscViewerFormat format; 561e6985dafSRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 56249bd79ccSDebojyoti Ghosh Mat B; 563e6985dafSRichard Tran Mills PetscInt i; 564e6985dafSRichard Tran Mills PetscBool ismpikaij; 56549bd79ccSDebojyoti Ghosh 56649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij)); 5685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetFormat(viewer,&format)); 569e6985dafSRichard Tran Mills if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) { 5705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"S and T have %" PetscInt_FMT " rows and %" PetscInt_FMT " columns\n",a->p,a->q)); 571e6985dafSRichard Tran Mills 572e6985dafSRichard Tran Mills /* Print appropriate details for S. */ 573e6985dafSRichard Tran Mills if (!a->S) { 5745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"S is NULL\n")); 575e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 5765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Entries of S are ")); 577e6985dafSRichard Tran Mills for (i=0; i<(a->p * a->q); i++) { 578e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 5795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->S[i]),(double)PetscImaginaryPart(a->S[i]))); 580e6985dafSRichard Tran Mills #else 5815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->S[i]))); 582e6985dafSRichard Tran Mills #endif 583e6985dafSRichard Tran Mills } 5845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n")); 58549bd79ccSDebojyoti Ghosh } 58649bd79ccSDebojyoti Ghosh 587e6985dafSRichard Tran Mills /* Print appropriate details for T. */ 588e6985dafSRichard Tran Mills if (a->isTI) { 5895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"T is the identity matrix\n")); 590e6985dafSRichard Tran Mills } else if (!a->T) { 5915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"T is NULL\n")); 592e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 5935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Entries of T are ")); 594e6985dafSRichard Tran Mills for (i=0; i<(a->p * a->q); i++) { 595e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 5965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->T[i]),(double)PetscImaginaryPart(a->T[i]))); 597e6985dafSRichard Tran Mills #else 5985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->T[i]))); 599e6985dafSRichard Tran Mills #endif 600e6985dafSRichard Tran Mills } 6015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"\n")); 602e6985dafSRichard Tran Mills } 60349bd79ccSDebojyoti Ghosh 604e6985dafSRichard Tran Mills /* Now print details for the AIJ matrix, using the AIJ viewer. */ 6055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Now viewing the associated AIJ matrix:\n")); 606e6985dafSRichard Tran Mills if (ismpikaij) { 607e6985dafSRichard Tran Mills Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 6085f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(b->A,viewer)); 609e6985dafSRichard Tran Mills } else { 6105f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(a->AIJ,viewer)); 611e6985dafSRichard Tran Mills } 612e6985dafSRichard Tran Mills 613e6985dafSRichard Tran Mills } else { 614e6985dafSRichard Tran Mills /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */ 6155f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 6165f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,viewer)); 6175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 618e6985dafSRichard Tran Mills } 61949bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 62049bd79ccSDebojyoti Ghosh } 62149bd79ccSDebojyoti Ghosh 62249bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_MPIKAIJ(Mat A) 62349bd79ccSDebojyoti Ghosh { 62449bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 62549bd79ccSDebojyoti Ghosh 62649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 6275f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&b->AIJ)); 6285f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&b->OAIJ)); 6295f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&b->A)); 6305f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&b->ctx)); 6315f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b->w)); 6325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->S)); 6335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->T)); 6345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b->ibdiag)); 6355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatGetDiagonalBlock_C",NULL)); 6365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatConvert_mpikaij_mpiaij_C",NULL)); 6375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(A->data)); 63849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 63949bd79ccSDebojyoti Ghosh } 64049bd79ccSDebojyoti Ghosh 64149bd79ccSDebojyoti Ghosh /* --------------------------------------------------------------------------------------*/ 64249bd79ccSDebojyoti Ghosh 64349bd79ccSDebojyoti Ghosh /* zz = yy + Axx */ 644836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_SeqKAIJ(Mat A,Vec xx,Vec yy,Vec zz) 64549bd79ccSDebojyoti Ghosh { 64649bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 64749bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 64849bd79ccSDebojyoti Ghosh const PetscScalar *s = b->S, *t = b->T; 64949bd79ccSDebojyoti Ghosh const PetscScalar *x,*v,*bx; 65049bd79ccSDebojyoti Ghosh PetscScalar *y,*sums; 65149bd79ccSDebojyoti Ghosh const PetscInt m = b->AIJ->rmap->n,*idx,*ii; 65249bd79ccSDebojyoti Ghosh PetscInt n,i,jrow,j,l,p=b->p,q=b->q,k; 65349bd79ccSDebojyoti Ghosh 65449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 65549bd79ccSDebojyoti Ghosh if (!yy) { 6565f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(zz,0.0)); 65749bd79ccSDebojyoti Ghosh } else { 6585f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(yy,zz)); 65949bd79ccSDebojyoti Ghosh } 66049bd79ccSDebojyoti Ghosh if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(0); 66149bd79ccSDebojyoti Ghosh 6625f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(xx,&x)); 6635f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(zz,&y)); 66449bd79ccSDebojyoti Ghosh idx = a->j; 66549bd79ccSDebojyoti Ghosh v = a->a; 66649bd79ccSDebojyoti Ghosh ii = a->i; 66749bd79ccSDebojyoti Ghosh 66849bd79ccSDebojyoti Ghosh if (b->isTI) { 66949bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 67049bd79ccSDebojyoti Ghosh jrow = ii[i]; 67149bd79ccSDebojyoti Ghosh n = ii[i+1] - jrow; 67249bd79ccSDebojyoti Ghosh sums = y + p*i; 67349bd79ccSDebojyoti Ghosh for (j=0; j<n; j++) { 67449bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 67549bd79ccSDebojyoti Ghosh sums[k] += v[jrow+j]*x[q*idx[jrow+j]+k]; 67649bd79ccSDebojyoti Ghosh } 67749bd79ccSDebojyoti Ghosh } 67849bd79ccSDebojyoti Ghosh } 6795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(3.0*(a->nz)*p)); 68049bd79ccSDebojyoti Ghosh } else if (t) { 68149bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 68249bd79ccSDebojyoti Ghosh jrow = ii[i]; 68349bd79ccSDebojyoti Ghosh n = ii[i+1] - jrow; 68449bd79ccSDebojyoti Ghosh sums = y + p*i; 68549bd79ccSDebojyoti Ghosh for (j=0; j<n; j++) { 68649bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 68749bd79ccSDebojyoti Ghosh for (l=0; l<q; l++) { 68849bd79ccSDebojyoti Ghosh sums[k] += v[jrow+j]*t[k+l*p]*x[q*idx[jrow+j]+l]; 68949bd79ccSDebojyoti Ghosh } 69049bd79ccSDebojyoti Ghosh } 69149bd79ccSDebojyoti Ghosh } 69249bd79ccSDebojyoti Ghosh } 693234d9204SRichard Tran Mills /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do), 694234d9204SRichard Tran Mills * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T), 695234d9204SRichard Tran Mills * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter 696234d9204SRichard Tran Mills * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */ 6975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops((2.0*p*q-p)*m+2.0*p*a->nz)); 69849bd79ccSDebojyoti Ghosh } 69949bd79ccSDebojyoti Ghosh if (s) { 70049bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 70149bd79ccSDebojyoti Ghosh sums = y + p*i; 70249bd79ccSDebojyoti Ghosh bx = x + q*i; 70349bd79ccSDebojyoti Ghosh if (i < b->AIJ->cmap->n) { 70449bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 70549bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 70649bd79ccSDebojyoti Ghosh sums[k] += s[k+j*p]*bx[j]; 70749bd79ccSDebojyoti Ghosh } 70849bd79ccSDebojyoti Ghosh } 70949bd79ccSDebojyoti Ghosh } 71049bd79ccSDebojyoti Ghosh } 7115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(2.0*m*p*q)); 71249bd79ccSDebojyoti Ghosh } 71349bd79ccSDebojyoti Ghosh 7145f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(xx,&x)); 7155f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(zz,&y)); 71649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 71749bd79ccSDebojyoti Ghosh } 71849bd79ccSDebojyoti Ghosh 719bb6fb833SRichard Tran Mills PetscErrorCode MatMult_SeqKAIJ(Mat A,Vec xx,Vec yy) 72049bd79ccSDebojyoti Ghosh { 72149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 7225f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd_SeqKAIJ(A,xx,PETSC_NULL,yy)); 72349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 72449bd79ccSDebojyoti Ghosh } 72549bd79ccSDebojyoti Ghosh 72649bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h> 72749bd79ccSDebojyoti Ghosh 728bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A,const PetscScalar **values) 72949bd79ccSDebojyoti Ghosh { 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; 7402c71b3e2SJacob Faibussowitsch PetscCheckFalse(p != q,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATKAIJ: Block size must be square to calculate inverse."); 7412c71b3e2SJacob Faibussowitsch PetscCheckFalse((!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; 74849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 74949bd79ccSDebojyoti Ghosh } 75049bd79ccSDebojyoti Ghosh if (!b->ibdiag) { 7515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(dof2*m,&b->ibdiag)); 7525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectMemory((PetscObject)A,dof2*m*sizeof(PetscScalar))); 75349bd79ccSDebojyoti Ghosh } 75449bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 75549bd79ccSDebojyoti Ghosh diag = b->ibdiag; 75649bd79ccSDebojyoti Ghosh 7575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(dof,&v_work,dof,&v_pivots)); 75849bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 75949bd79ccSDebojyoti Ghosh if (S) { 7605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(diag,S,dof2*sizeof(PetscScalar))); 76149bd79ccSDebojyoti Ghosh } else { 7625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(diag,dof2*sizeof(PetscScalar))); 76349bd79ccSDebojyoti Ghosh } 76449bd79ccSDebojyoti Ghosh if (b->isTI) { 76549bd79ccSDebojyoti Ghosh aval = 0; 76649bd79ccSDebojyoti Ghosh for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j]; 76749bd79ccSDebojyoti Ghosh for (j=0; j<dof; j++) diag[j+dof*j] += aval; 76849bd79ccSDebojyoti Ghosh } else if (T) { 76949bd79ccSDebojyoti Ghosh aval = 0; 77049bd79ccSDebojyoti Ghosh for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j]; 77149bd79ccSDebojyoti Ghosh for (j=0; j<dof2; j++) diag[j] += aval*T[j]; 77249bd79ccSDebojyoti Ghosh } 7735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscKernel_A_gets_inverse_A(dof,diag,v_pivots,v_work,PETSC_FALSE,NULL)); 77449bd79ccSDebojyoti Ghosh diag += dof2; 77549bd79ccSDebojyoti Ghosh } 7765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(v_work,v_pivots)); 77749bd79ccSDebojyoti Ghosh 77849bd79ccSDebojyoti Ghosh b->ibdiagvalid = PETSC_TRUE; 77949bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 78049bd79ccSDebojyoti Ghosh } 78149bd79ccSDebojyoti Ghosh 78249bd79ccSDebojyoti Ghosh static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A,Mat *B) 78349bd79ccSDebojyoti Ghosh { 78449bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ*) A->data; 78549bd79ccSDebojyoti Ghosh 78649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 78749bd79ccSDebojyoti Ghosh *B = kaij->AIJ; 78849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 78949bd79ccSDebojyoti Ghosh } 79049bd79ccSDebojyoti Ghosh 791fac53fa1SPierre Jolivet static PetscErrorCode MatConvert_KAIJ_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 792fac53fa1SPierre Jolivet { 793fac53fa1SPierre Jolivet Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 794fac53fa1SPierre Jolivet Mat AIJ,OAIJ,B; 795fac53fa1SPierre Jolivet PetscInt *d_nnz,*o_nnz = NULL,nz,i,j,m,d; 796fac53fa1SPierre Jolivet const PetscInt p = a->p,q = a->q; 797fac53fa1SPierre Jolivet PetscBool ismpikaij,missing; 798fac53fa1SPierre Jolivet 799fac53fa1SPierre Jolivet PetscFunctionBegin; 800fac53fa1SPierre Jolivet if (reuse != MAT_REUSE_MATRIX) { 8015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij)); 802fac53fa1SPierre Jolivet if (ismpikaij) { 803fac53fa1SPierre Jolivet Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 804fac53fa1SPierre Jolivet AIJ = ((Mat_SeqKAIJ*)b->AIJ->data)->AIJ; 805fac53fa1SPierre Jolivet OAIJ = ((Mat_SeqKAIJ*)b->OAIJ->data)->AIJ; 806fac53fa1SPierre Jolivet } else { 807fac53fa1SPierre Jolivet AIJ = a->AIJ; 808fac53fa1SPierre Jolivet OAIJ = NULL; 809fac53fa1SPierre Jolivet } 8105f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&B)); 8115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 8125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATAIJ)); 8135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(AIJ,&m,NULL)); 8145f80ce2aSJacob Faibussowitsch CHKERRQ(MatMissingDiagonal(AIJ,&missing,&d)); /* assumption that all successive rows will have a missing diagonal */ 815fac53fa1SPierre Jolivet if (!missing || !a->S) d = m; 8165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*p,&d_nnz)); 817fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8185f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow_SeqAIJ(AIJ,i,&nz,NULL,NULL)); 819fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) d_nnz[i*p + j] = nz*q + (i >= d)*q; 8205f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow_SeqAIJ(AIJ,i,&nz,NULL,NULL)); 821fac53fa1SPierre Jolivet } 822fac53fa1SPierre Jolivet if (OAIJ) { 8235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m*p,&o_nnz)); 824fac53fa1SPierre Jolivet for (i = 0; i < m; ++i) { 8255f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow_SeqAIJ(OAIJ,i,&nz,NULL,NULL)); 826fac53fa1SPierre Jolivet for (j = 0; j < p; ++j) o_nnz[i*p + j] = nz*q; 8275f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow_SeqAIJ(OAIJ,i,&nz,NULL,NULL)); 828fac53fa1SPierre Jolivet } 8295f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(B,0,d_nnz,0,o_nnz)); 830fac53fa1SPierre Jolivet } else { 8315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(B,0,d_nnz)); 832fac53fa1SPierre Jolivet } 8335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(d_nnz)); 8345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(o_nnz)); 835fac53fa1SPierre Jolivet } else B = *newmat; 8365f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert_Basic(A,newtype,MAT_REUSE_MATRIX,&B)); 837fac53fa1SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 8385f80ce2aSJacob Faibussowitsch CHKERRQ(MatHeaderReplace(A,&B)); 839fac53fa1SPierre Jolivet } else *newmat = B; 840fac53fa1SPierre Jolivet PetscFunctionReturn(0); 841fac53fa1SPierre Jolivet } 842fac53fa1SPierre Jolivet 84349bd79ccSDebojyoti Ghosh PetscErrorCode MatSOR_SeqKAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 84449bd79ccSDebojyoti Ghosh { 84549bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 84649bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *kaij = (Mat_SeqKAIJ*) A->data; 84749bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ*)kaij->AIJ->data; 84849bd79ccSDebojyoti Ghosh const PetscScalar *aa = a->a, *T = kaij->T, *v; 84949bd79ccSDebojyoti Ghosh const PetscInt m = kaij->AIJ->rmap->n, *ai=a->i, *aj=a->j, p = kaij->p, q = kaij->q, *diag, *vi; 85049bd79ccSDebojyoti Ghosh const PetscScalar *b, *xb, *idiag; 85149bd79ccSDebojyoti Ghosh PetscScalar *x, *work, *workt, *w, *y, *arr, *t, *arrt; 85249bd79ccSDebojyoti Ghosh PetscInt i, j, k, i2, bs, bs2, nz; 85349bd79ccSDebojyoti Ghosh 85449bd79ccSDebojyoti Ghosh PetscFunctionBegin; 85549bd79ccSDebojyoti Ghosh its = its*lits; 8562c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag & SOR_EISENSTAT,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 8572c71b3e2SJacob Faibussowitsch PetscCheckFalse(its <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 858*28b400f6SJacob Faibussowitsch PetscCheck(!fshift,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for diagonal shift"); 8592c71b3e2SJacob Faibussowitsch PetscCheckFalse((flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER),PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for applying upper or lower triangular parts"); 8602c71b3e2SJacob Faibussowitsch PetscCheckFalse(p != q,PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSOR for KAIJ: No support for non-square dense blocks"); 86149bd79ccSDebojyoti Ghosh else {bs = p; bs2 = bs*bs; } 86249bd79ccSDebojyoti Ghosh 86349bd79ccSDebojyoti Ghosh if (!m) PetscFunctionReturn(0); 86449bd79ccSDebojyoti Ghosh 8655f80ce2aSJacob Faibussowitsch if (!kaij->ibdiagvalid) CHKERRQ(MatInvertBlockDiagonal_SeqKAIJ(A,NULL)); 86649bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 86749bd79ccSDebojyoti Ghosh diag = a->diag; 86849bd79ccSDebojyoti Ghosh 86949bd79ccSDebojyoti Ghosh if (!kaij->sor.setup) { 8705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc5(bs,&kaij->sor.w,bs,&kaij->sor.y,m*bs,&kaij->sor.work,m*bs,&kaij->sor.t,m*bs2,&kaij->sor.arr)); 87149bd79ccSDebojyoti Ghosh kaij->sor.setup = PETSC_TRUE; 87249bd79ccSDebojyoti Ghosh } 87349bd79ccSDebojyoti Ghosh y = kaij->sor.y; 87449bd79ccSDebojyoti Ghosh w = kaij->sor.w; 87549bd79ccSDebojyoti Ghosh work = kaij->sor.work; 87649bd79ccSDebojyoti Ghosh t = kaij->sor.t; 87749bd79ccSDebojyoti Ghosh arr = kaij->sor.arr; 87849bd79ccSDebojyoti Ghosh 87949bd79ccSDebojyoti Ghosh ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 8805f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(bb,&b)); 88149bd79ccSDebojyoti Ghosh 88249bd79ccSDebojyoti Ghosh if (flag & SOR_ZERO_INITIAL_GUESS) { 88349bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 88449bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,b,idiag,x); /* x[0:bs] <- D^{-1} b[0:bs] */ 8855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(t,b,bs*sizeof(PetscScalar))); 88649bd79ccSDebojyoti Ghosh i2 = bs; 88749bd79ccSDebojyoti Ghosh idiag += bs2; 88849bd79ccSDebojyoti Ghosh for (i=1; i<m; i++) { 88949bd79ccSDebojyoti Ghosh v = aa + ai[i]; 89049bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 89149bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 89249bd79ccSDebojyoti Ghosh 89349bd79ccSDebojyoti Ghosh if (T) { /* b - T (Arow * x) */ 8945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(w,bs*sizeof(PetscScalar))); 89549bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 89649bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k]; 89749bd79ccSDebojyoti Ghosh } 89849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs,w,T,&t[i2]); 89949bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) t[i2+k] += b[i2+k]; 90049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar))); 90249bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 90349bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) t[i2+k] -= v[j] * x[vi[j]*bs+k]; 90449bd79ccSDebojyoti Ghosh } 90549bd79ccSDebojyoti Ghosh } else { 9065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar))); 90749bd79ccSDebojyoti Ghosh } 90849bd79ccSDebojyoti Ghosh 90949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,t+i2,idiag,y); 91049bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) x[i2+j] = omega * y[j]; 91149bd79ccSDebojyoti Ghosh 91249bd79ccSDebojyoti Ghosh idiag += bs2; 91349bd79ccSDebojyoti Ghosh i2 += bs; 91449bd79ccSDebojyoti Ghosh } 91549bd79ccSDebojyoti Ghosh /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 9165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(1.0*bs2*a->nz)); 91749bd79ccSDebojyoti Ghosh xb = t; 91849bd79ccSDebojyoti Ghosh } else xb = b; 91949bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 92049bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag+bs2*(m-1); 92149bd79ccSDebojyoti Ghosh i2 = bs * (m-1); 9225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar))); 92349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2); 92449bd79ccSDebojyoti Ghosh i2 -= bs; 92549bd79ccSDebojyoti Ghosh idiag -= bs2; 92649bd79ccSDebojyoti Ghosh for (i=m-2; i>=0; i--) { 92749bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1 ; 92849bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 92949bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 93049bd79ccSDebojyoti Ghosh 93149bd79ccSDebojyoti Ghosh if (T) { /* FIXME: This branch untested */ 9325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar))); 93349bd79ccSDebojyoti Ghosh /* copy all rows of x that are needed into contiguous space */ 93449bd79ccSDebojyoti Ghosh workt = work; 93549bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 9365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 93749bd79ccSDebojyoti Ghosh workt += bs; 93849bd79ccSDebojyoti Ghosh } 93949bd79ccSDebojyoti Ghosh arrt = arr; 94049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 9415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 94249bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 94349bd79ccSDebojyoti Ghosh arrt += bs2; 94449bd79ccSDebojyoti Ghosh } 94549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 94649bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,t+i2,bs*sizeof(PetscScalar))); 94849bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 94949bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k]; 95049bd79ccSDebojyoti Ghosh } 95149bd79ccSDebojyoti Ghosh } 95249bd79ccSDebojyoti Ghosh 95349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); /* RHS incorrect for omega != 1.0 */ 95449bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) x[i2+j] = (1.0-omega) * x[i2+j] + omega * y[j]; 95549bd79ccSDebojyoti Ghosh 95649bd79ccSDebojyoti Ghosh idiag -= bs2; 95749bd79ccSDebojyoti Ghosh i2 -= bs; 95849bd79ccSDebojyoti Ghosh } 9595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(1.0*bs2*(a->nz))); 96049bd79ccSDebojyoti Ghosh } 96149bd79ccSDebojyoti Ghosh its--; 96249bd79ccSDebojyoti Ghosh } 96349bd79ccSDebojyoti Ghosh while (its--) { /* FIXME: This branch not updated */ 96449bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 96549bd79ccSDebojyoti Ghosh i2 = 0; 96649bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 96749bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 9685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar))); 96949bd79ccSDebojyoti Ghosh 97049bd79ccSDebojyoti Ghosh v = aa + ai[i]; 97149bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 97249bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 97349bd79ccSDebojyoti Ghosh workt = work; 97449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 9755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 97649bd79ccSDebojyoti Ghosh workt += bs; 97749bd79ccSDebojyoti Ghosh } 97849bd79ccSDebojyoti Ghosh arrt = arr; 97949bd79ccSDebojyoti Ghosh if (T) { 98049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 9815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 98249bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 98349bd79ccSDebojyoti Ghosh arrt += bs2; 98449bd79ccSDebojyoti Ghosh } 98549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 98649bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 98749bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 9885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(arrt,bs2*sizeof(PetscScalar))); 98949bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 99049bd79ccSDebojyoti Ghosh arrt += bs2; 99149bd79ccSDebojyoti Ghosh } 99249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 99349bd79ccSDebojyoti Ghosh } 9945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(t+i2,w,bs*sizeof(PetscScalar))); 99549bd79ccSDebojyoti Ghosh 99649bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 99749bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 99849bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 99949bd79ccSDebojyoti Ghosh workt = work; 100049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 100249bd79ccSDebojyoti Ghosh workt += bs; 100349bd79ccSDebojyoti Ghosh } 100449bd79ccSDebojyoti Ghosh arrt = arr; 100549bd79ccSDebojyoti Ghosh if (T) { 100649bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 100849bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 100949bd79ccSDebojyoti Ghosh arrt += bs2; 101049bd79ccSDebojyoti Ghosh } 101149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 101249bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 101349bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(arrt,bs2*sizeof(PetscScalar))); 101549bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 101649bd79ccSDebojyoti Ghosh arrt += bs2; 101749bd79ccSDebojyoti Ghosh } 101849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 101949bd79ccSDebojyoti Ghosh } 102049bd79ccSDebojyoti Ghosh 102149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 102249bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 102349bd79ccSDebojyoti Ghosh 102449bd79ccSDebojyoti Ghosh idiag += bs2; 102549bd79ccSDebojyoti Ghosh i2 += bs; 102649bd79ccSDebojyoti Ghosh } 102749bd79ccSDebojyoti Ghosh xb = t; 102849bd79ccSDebojyoti Ghosh } 102949bd79ccSDebojyoti Ghosh else xb = b; 103049bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 103149bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag+bs2*(m-1); 103249bd79ccSDebojyoti Ghosh i2 = bs * (m-1); 103349bd79ccSDebojyoti Ghosh if (xb == b) { 103449bd79ccSDebojyoti Ghosh for (i=m-1; i>=0; i--) { 10355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar))); 103649bd79ccSDebojyoti Ghosh 103749bd79ccSDebojyoti Ghosh v = aa + ai[i]; 103849bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 103949bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 104049bd79ccSDebojyoti Ghosh workt = work; 104149bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 104349bd79ccSDebojyoti Ghosh workt += bs; 104449bd79ccSDebojyoti Ghosh } 104549bd79ccSDebojyoti Ghosh arrt = arr; 104649bd79ccSDebojyoti Ghosh if (T) { 104749bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 104949bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 105049bd79ccSDebojyoti Ghosh arrt += bs2; 105149bd79ccSDebojyoti Ghosh } 105249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 105349bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 105449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(arrt,bs2*sizeof(PetscScalar))); 105649bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*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 } 106149bd79ccSDebojyoti Ghosh 106249bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 106349bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 106449bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 106549bd79ccSDebojyoti Ghosh workt = work; 106649bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 106849bd79ccSDebojyoti Ghosh workt += bs; 106949bd79ccSDebojyoti Ghosh } 107049bd79ccSDebojyoti Ghosh arrt = arr; 107149bd79ccSDebojyoti Ghosh if (T) { 107249bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 107449bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 107549bd79ccSDebojyoti Ghosh arrt += bs2; 107649bd79ccSDebojyoti Ghosh } 107749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 107849bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 107949bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(arrt,bs2*sizeof(PetscScalar))); 108149bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*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 } 108649bd79ccSDebojyoti Ghosh 108749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 108849bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 108949bd79ccSDebojyoti Ghosh } 109049bd79ccSDebojyoti Ghosh } else { 109149bd79ccSDebojyoti Ghosh for (i=m-1; i>=0; i--) { 10925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar))); 109349bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 109449bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 109549bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 109649bd79ccSDebojyoti Ghosh workt = work; 109749bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 10985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar))); 109949bd79ccSDebojyoti Ghosh workt += bs; 110049bd79ccSDebojyoti Ghosh } 110149bd79ccSDebojyoti Ghosh arrt = arr; 110249bd79ccSDebojyoti Ghosh if (T) { 110349bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 11045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar))); 110549bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 110649bd79ccSDebojyoti Ghosh arrt += bs2; 110749bd79ccSDebojyoti Ghosh } 110849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 110949bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 111049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 11115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(arrt,bs2*sizeof(PetscScalar))); 111249bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*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 } 111749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 111849bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 111949bd79ccSDebojyoti Ghosh } 112049bd79ccSDebojyoti Ghosh } 11215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(1.0*bs2*(a->nz))); 112249bd79ccSDebojyoti Ghosh } 112349bd79ccSDebojyoti Ghosh } 112449bd79ccSDebojyoti Ghosh 112549bd79ccSDebojyoti Ghosh ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 11265f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(bb,&b)); 112749bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 112849bd79ccSDebojyoti Ghosh } 112949bd79ccSDebojyoti Ghosh 113049bd79ccSDebojyoti Ghosh /*===================================================================================*/ 113149bd79ccSDebojyoti Ghosh 1132836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_MPIKAIJ(Mat A,Vec xx,Vec yy,Vec zz) 113349bd79ccSDebojyoti Ghosh { 113449bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 113549bd79ccSDebojyoti Ghosh 113649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 113749bd79ccSDebojyoti Ghosh if (!yy) { 11385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(zz,0.0)); 113949bd79ccSDebojyoti Ghosh } else { 11405f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(yy,zz)); 114149bd79ccSDebojyoti Ghosh } 11425f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */ 114349bd79ccSDebojyoti Ghosh /* start the scatter */ 11445f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD)); 11455f80ce2aSJacob Faibussowitsch CHKERRQ((*b->AIJ->ops->multadd)(b->AIJ,xx,zz,zz)); 11465f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD)); 11475f80ce2aSJacob Faibussowitsch CHKERRQ((*b->OAIJ->ops->multadd)(b->OAIJ,b->w,zz,zz)); 114849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 114949bd79ccSDebojyoti Ghosh } 115049bd79ccSDebojyoti Ghosh 1151bb6fb833SRichard Tran Mills PetscErrorCode MatMult_MPIKAIJ(Mat A,Vec xx,Vec yy) 115249bd79ccSDebojyoti Ghosh { 115349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 11545f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd_MPIKAIJ(A,xx,PETSC_NULL,yy)); 115549bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 115649bd79ccSDebojyoti Ghosh } 115749bd79ccSDebojyoti Ghosh 1158bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A,const PetscScalar **values) 115949bd79ccSDebojyoti Ghosh { 116049bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 116149bd79ccSDebojyoti Ghosh 116249bd79ccSDebojyoti Ghosh PetscFunctionBegin; 11635f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ is up to date. */ 11645f80ce2aSJacob Faibussowitsch CHKERRQ((*b->AIJ->ops->invertblockdiagonal)(b->AIJ,values)); 116549bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 116649bd79ccSDebojyoti Ghosh } 116749bd79ccSDebojyoti Ghosh 116849bd79ccSDebojyoti Ghosh /* ----------------------------------------------------------------*/ 116949bd79ccSDebojyoti Ghosh 117049bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values) 117149bd79ccSDebojyoti Ghosh { 117249bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*) A->data; 11731ca5ffdbSRichard Tran Mills PetscErrorCode diag = PETSC_FALSE; 117449bd79ccSDebojyoti Ghosh PetscInt nzaij,nz,*colsaij,*idx,i,j,p=b->p,q=b->q,r=row/p,s=row%p,c; 117549bd79ccSDebojyoti Ghosh PetscScalar *vaij,*v,*S=b->S,*T=b->T; 117649bd79ccSDebojyoti Ghosh 117749bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1178*28b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active"); 117949bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 11802c71b3e2SJacob Faibussowitsch PetscCheckFalse(row < 0 || row >= A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " out of range",row); 118149bd79ccSDebojyoti Ghosh 118249bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 118349bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 118449bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 118549bd79ccSDebojyoti Ghosh if (values) *values = NULL; 118649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 118749bd79ccSDebojyoti Ghosh } 118849bd79ccSDebojyoti Ghosh 118949bd79ccSDebojyoti Ghosh if (T || b->isTI) { 11905f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow_SeqAIJ(b->AIJ,r,&nzaij,&colsaij,&vaij)); 119149bd79ccSDebojyoti Ghosh c = nzaij; 119249bd79ccSDebojyoti Ghosh for (i=0; i<nzaij; i++) { 119349bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 119449bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 119549bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 119649bd79ccSDebojyoti Ghosh c = i; 119749bd79ccSDebojyoti Ghosh } 119849bd79ccSDebojyoti Ghosh } 119949bd79ccSDebojyoti Ghosh } else nzaij = c = 0; 120049bd79ccSDebojyoti Ghosh 120149bd79ccSDebojyoti Ghosh /* calculate size of row */ 120249bd79ccSDebojyoti Ghosh nz = 0; 120349bd79ccSDebojyoti Ghosh if (S) nz += q; 120449bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (nzaij-1)*q : nzaij*q); 120549bd79ccSDebojyoti Ghosh 120649bd79ccSDebojyoti Ghosh if (cols || values) { 12075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nz,&idx,nz,&v)); 120838822f9dSRichard Tran Mills for (i=0; i<q; i++) { 120938822f9dSRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 121038822f9dSRichard Tran Mills v[i] = 0.0; 121138822f9dSRichard Tran Mills } 121249bd79ccSDebojyoti Ghosh if (b->isTI) { 121349bd79ccSDebojyoti Ghosh for (i=0; i<nzaij; i++) { 121449bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 121549bd79ccSDebojyoti Ghosh idx[i*q+j] = colsaij[i]*q+j; 121649bd79ccSDebojyoti Ghosh v[i*q+j] = (j==s ? vaij[i] : 0); 121749bd79ccSDebojyoti Ghosh } 121849bd79ccSDebojyoti Ghosh } 121949bd79ccSDebojyoti Ghosh } else if (T) { 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] = vaij[i]*T[s+j*p]; 122449bd79ccSDebojyoti Ghosh } 122549bd79ccSDebojyoti Ghosh } 122649bd79ccSDebojyoti Ghosh } 122749bd79ccSDebojyoti Ghosh if (S) { 122849bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 122949bd79ccSDebojyoti Ghosh idx[c*q+j] = r*q+j; 123049bd79ccSDebojyoti Ghosh v[c*q+j] += S[s+j*p]; 123149bd79ccSDebojyoti Ghosh } 123249bd79ccSDebojyoti Ghosh } 123349bd79ccSDebojyoti Ghosh } 123449bd79ccSDebojyoti Ghosh 123549bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 123649bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 123749bd79ccSDebojyoti Ghosh if (values) *values = v; 123849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 123949bd79ccSDebojyoti Ghosh } 124049bd79ccSDebojyoti Ghosh 124149bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 124249bd79ccSDebojyoti Ghosh { 124349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1244cb4a9cd9SHong Zhang if (nz) *nz = 0; 12455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(*idx,*v)); 124649bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE; 124749bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 124849bd79ccSDebojyoti Ghosh } 124949bd79ccSDebojyoti Ghosh 125049bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values) 125149bd79ccSDebojyoti Ghosh { 125249bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*) A->data; 125349bd79ccSDebojyoti Ghosh Mat AIJ = b->A; 1254fc64b2cfSRichard Tran Mills PetscBool diag = PETSC_FALSE; 1255761d359dSRichard Tran Mills Mat MatAIJ,MatOAIJ; 125649bd79ccSDebojyoti Ghosh const PetscInt rstart=A->rmap->rstart,rend=A->rmap->rend,p=b->p,q=b->q,*garray; 1257389eba51SJed Brown PetscInt nz,*idx,ncolsaij = 0,ncolsoaij = 0,*colsaij,*colsoaij,r,s,c,i,j,lrow; 125849bd79ccSDebojyoti Ghosh PetscScalar *v,*vals,*ovals,*S=b->S,*T=b->T; 125949bd79ccSDebojyoti Ghosh 126049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 12615f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJ_build_AIJ_OAIJ(A)); /* Ensure b->AIJ and b->OAIJ are up to date. */ 1262761d359dSRichard Tran Mills MatAIJ = ((Mat_SeqKAIJ*)b->AIJ->data)->AIJ; 1263761d359dSRichard Tran Mills MatOAIJ = ((Mat_SeqKAIJ*)b->OAIJ->data)->AIJ; 1264*28b400f6SJacob Faibussowitsch PetscCheck(!b->getrowactive,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active"); 126549bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 12662c71b3e2SJacob Faibussowitsch PetscCheckFalse(row < rstart || row >= rend,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only local rows"); 126749bd79ccSDebojyoti Ghosh lrow = row - rstart; 126849bd79ccSDebojyoti Ghosh 126949bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 127049bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 127149bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 127249bd79ccSDebojyoti Ghosh if (values) *values = NULL; 127349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 127449bd79ccSDebojyoti Ghosh } 127549bd79ccSDebojyoti Ghosh 127649bd79ccSDebojyoti Ghosh r = lrow/p; 127749bd79ccSDebojyoti Ghosh s = lrow%p; 127849bd79ccSDebojyoti Ghosh 127949bd79ccSDebojyoti Ghosh if (T || b->isTI) { 12805f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(AIJ,NULL,NULL,&garray)); 12815f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow_SeqAIJ(MatAIJ,lrow/p,&ncolsaij,&colsaij,&vals)); 12825f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow_SeqAIJ(MatOAIJ,lrow/p,&ncolsoaij,&colsoaij,&ovals)); 128349bd79ccSDebojyoti Ghosh 128449bd79ccSDebojyoti Ghosh c = ncolsaij + ncolsoaij; 128549bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 128649bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 128749bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 128849bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 128949bd79ccSDebojyoti Ghosh c = i; 129049bd79ccSDebojyoti Ghosh } 129149bd79ccSDebojyoti Ghosh } 129249bd79ccSDebojyoti Ghosh } else c = 0; 129349bd79ccSDebojyoti Ghosh 129449bd79ccSDebojyoti Ghosh /* calculate size of row */ 129549bd79ccSDebojyoti Ghosh nz = 0; 129649bd79ccSDebojyoti Ghosh if (S) nz += q; 129749bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (ncolsaij+ncolsoaij-1)*q : (ncolsaij+ncolsoaij)*q); 129849bd79ccSDebojyoti Ghosh 129949bd79ccSDebojyoti Ghosh if (cols || values) { 13005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(nz,&idx,nz,&v)); 1301a437a796SRichard Tran Mills for (i=0; i<q; i++) { 1302a437a796SRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 1303a437a796SRichard Tran Mills v[i] = 0.0; 1304a437a796SRichard Tran Mills } 130549bd79ccSDebojyoti Ghosh if (b->isTI) { 130649bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 130749bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 130849bd79ccSDebojyoti Ghosh idx[i*q+j] = (colsaij[i]+rstart/p)*q+j; 130949bd79ccSDebojyoti Ghosh v[i*q+j] = (j==s ? vals[i] : 0.0); 131049bd79ccSDebojyoti Ghosh } 131149bd79ccSDebojyoti Ghosh } 131249bd79ccSDebojyoti Ghosh for (i=0; i<ncolsoaij; i++) { 131349bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 131449bd79ccSDebojyoti Ghosh idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j; 131549bd79ccSDebojyoti Ghosh v[(i+ncolsaij)*q+j] = (j==s ? ovals[i]: 0.0); 131649bd79ccSDebojyoti Ghosh } 131749bd79ccSDebojyoti Ghosh } 131849bd79ccSDebojyoti Ghosh } else if (T) { 131949bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 132049bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 132149bd79ccSDebojyoti Ghosh idx[i*q+j] = (colsaij[i]+rstart/p)*q+j; 132249bd79ccSDebojyoti Ghosh v[i*q+j] = vals[i]*T[s+j*p]; 132349bd79ccSDebojyoti Ghosh } 132449bd79ccSDebojyoti Ghosh } 132549bd79ccSDebojyoti Ghosh for (i=0; i<ncolsoaij; i++) { 132649bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 132749bd79ccSDebojyoti Ghosh idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j; 132849bd79ccSDebojyoti Ghosh v[(i+ncolsaij)*q+j] = ovals[i]*T[s+j*p]; 132949bd79ccSDebojyoti Ghosh } 133049bd79ccSDebojyoti Ghosh } 133149bd79ccSDebojyoti Ghosh } 133249bd79ccSDebojyoti Ghosh if (S) { 133349bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 133449bd79ccSDebojyoti Ghosh idx[c*q+j] = (r+rstart/p)*q+j; 133549bd79ccSDebojyoti Ghosh v[c*q+j] += S[s+j*p]; 133649bd79ccSDebojyoti Ghosh } 133749bd79ccSDebojyoti Ghosh } 133849bd79ccSDebojyoti Ghosh } 133949bd79ccSDebojyoti Ghosh 134049bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 134149bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 134249bd79ccSDebojyoti Ghosh if (values) *values = v; 134349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 134449bd79ccSDebojyoti Ghosh } 134549bd79ccSDebojyoti Ghosh 134649bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 134749bd79ccSDebojyoti Ghosh { 134849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 13495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(*idx,*v)); 135049bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE; 135149bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 135249bd79ccSDebojyoti Ghosh } 135349bd79ccSDebojyoti Ghosh 135449bd79ccSDebojyoti Ghosh PetscErrorCode MatCreateSubMatrix_KAIJ(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) 135549bd79ccSDebojyoti Ghosh { 135649bd79ccSDebojyoti Ghosh Mat A; 135749bd79ccSDebojyoti Ghosh 135849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 13595f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(mat,MATAIJ,MAT_INITIAL_MATRIX,&A)); 13605f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(A,isrow,iscol,cll,newmat)); 13615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 136249bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 136349bd79ccSDebojyoti Ghosh } 136449bd79ccSDebojyoti Ghosh 136549bd79ccSDebojyoti Ghosh /* ---------------------------------------------------------------------------------- */ 136649bd79ccSDebojyoti Ghosh /*@C 136749bd79ccSDebojyoti Ghosh MatCreateKAIJ - Creates a matrix type to be used for matrices of the following form: 136849bd79ccSDebojyoti Ghosh 136949bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 137049bd79ccSDebojyoti Ghosh 137149bd79ccSDebojyoti Ghosh where 137249bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 137349bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 137449bd79ccSDebojyoti Ghosh A is an AIJ (n \times n) matrix 137549bd79ccSDebojyoti Ghosh I is the identity matrix 137649bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 137749bd79ccSDebojyoti Ghosh 1378d3b90ce1SRichard Tran Mills S and T are always stored independently on all processes as PetscScalar arrays in column-major format. 137949bd79ccSDebojyoti Ghosh 138049bd79ccSDebojyoti Ghosh Collective 138149bd79ccSDebojyoti Ghosh 138249bd79ccSDebojyoti Ghosh Input Parameters: 138349bd79ccSDebojyoti Ghosh + A - the AIJ matrix 138449bd79ccSDebojyoti Ghosh . p - number of rows in S and T 1385d3b90ce1SRichard Tran Mills . q - number of columns in S and T 1386d3b90ce1SRichard Tran Mills . S - the S matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major) 1387d3b90ce1SRichard Tran Mills - T - the T matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major) 138849bd79ccSDebojyoti Ghosh 138949bd79ccSDebojyoti Ghosh Output Parameter: 139049bd79ccSDebojyoti Ghosh . kaij - the new KAIJ matrix 139149bd79ccSDebojyoti Ghosh 1392d3b90ce1SRichard Tran Mills Notes: 1393d3b90ce1SRichard Tran Mills This function increases the reference count on the AIJ matrix, so the user is free to destroy the matrix if it is not needed. 1394d3b90ce1SRichard Tran Mills Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix. 139549bd79ccSDebojyoti Ghosh 1396761d359dSRichard Tran Mills Developer Notes: 1397761d359dSRichard Tran Mills In the MATMPIKAIJ case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state 1398761d359dSRichard Tran Mills of the AIJ matrix 'A' that describes the blockwise action of the MATMPIKAIJ matrix and, if the object state has changed, lazily 1399761d359dSRichard Tran Mills rebuilding 'AIJ' and 'OAIJ' just before executing operations with the MATMPIKAIJ matrix. If new types of operations are added, 1400761d359dSRichard Tran Mills routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine). 1401761d359dSRichard Tran Mills 140249bd79ccSDebojyoti Ghosh Level: advanced 140349bd79ccSDebojyoti Ghosh 14040567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MATKAIJ 140549bd79ccSDebojyoti Ghosh @*/ 140649bd79ccSDebojyoti Ghosh PetscErrorCode MatCreateKAIJ(Mat A,PetscInt p,PetscInt q,const PetscScalar S[],const PetscScalar T[],Mat *kaij) 140749bd79ccSDebojyoti Ghosh { 140849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 14095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),kaij)); 14105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(*kaij,MATKAIJ)); 14115f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJSetAIJ(*kaij,A)); 14125f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJSetS(*kaij,p,q,S)); 14135f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJSetT(*kaij,p,q,T)); 14145f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(*kaij)); 14150567c835SRichard Tran Mills PetscFunctionReturn(0); 14160567c835SRichard Tran Mills } 141749bd79ccSDebojyoti Ghosh 14180567c835SRichard Tran Mills /*MC 14195881e567SRichard Tran Mills MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form 14205881e567SRichard Tran Mills [I \otimes S + A \otimes T], 14210567c835SRichard Tran Mills where 14225881e567SRichard Tran Mills S is a dense (p \times q) matrix, 14235881e567SRichard Tran Mills T is a dense (p \times q) matrix, 14245881e567SRichard Tran Mills A is an AIJ (n \times n) matrix, 14255881e567SRichard Tran Mills and I is the identity matrix. 14265881e567SRichard Tran Mills The resulting matrix is (np \times nq). 14270567c835SRichard Tran Mills 1428d3b90ce1SRichard Tran Mills S and T are always stored independently on all processes as PetscScalar arrays in column-major format. 14290567c835SRichard Tran Mills 14305881e567SRichard Tran Mills Notes: 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 14340567c835SRichard Tran Mills Level: advanced 14350567c835SRichard Tran Mills 14360567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MatCreateKAIJ() 14370567c835SRichard Tran Mills M*/ 14380567c835SRichard Tran Mills 14390567c835SRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A) 14400567c835SRichard Tran Mills { 14410567c835SRichard Tran Mills Mat_MPIKAIJ *b; 14420567c835SRichard Tran Mills PetscMPIInt size; 14430567c835SRichard Tran Mills 14440567c835SRichard Tran Mills PetscFunctionBegin; 14455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(A,&b)); 14460567c835SRichard Tran Mills A->data = (void*)b; 14470567c835SRichard Tran Mills 14485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemzero(A->ops,sizeof(struct _MatOps))); 14490567c835SRichard Tran Mills 1450f4259b30SLisandro Dalcin b->w = NULL; 14515f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 14520567c835SRichard Tran Mills if (size == 1) { 14535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)A,MATSEQKAIJ)); 14540567c835SRichard Tran Mills A->ops->destroy = MatDestroy_SeqKAIJ; 1455bb6fb833SRichard Tran Mills A->ops->mult = MatMult_SeqKAIJ; 1456bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_SeqKAIJ; 1457bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ; 14580567c835SRichard Tran Mills A->ops->getrow = MatGetRow_SeqKAIJ; 14590567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_SeqKAIJ; 14600567c835SRichard Tran Mills A->ops->sor = MatSOR_SeqKAIJ; 14615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqkaij_seqaij_C",MatConvert_KAIJ_AIJ)); 14620567c835SRichard Tran Mills } else { 14635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)A,MATMPIKAIJ)); 14640567c835SRichard Tran Mills A->ops->destroy = MatDestroy_MPIKAIJ; 1465bb6fb833SRichard Tran Mills A->ops->mult = MatMult_MPIKAIJ; 1466bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_MPIKAIJ; 1467bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ; 14680567c835SRichard Tran Mills A->ops->getrow = MatGetRow_MPIKAIJ; 14690567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_MPIKAIJ; 14705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIKAIJ)); 14715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatConvert_mpikaij_mpiaij_C",MatConvert_KAIJ_AIJ)); 14720567c835SRichard Tran Mills } 147306d61a37SPierre Jolivet A->ops->setup = MatSetUp_KAIJ; 147406d61a37SPierre Jolivet A->ops->view = MatView_KAIJ; 14750567c835SRichard Tran Mills A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ; 147649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 147749bd79ccSDebojyoti Ghosh } 1478