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 PetscErrorCode ierr; 5049bd79ccSDebojyoti Ghosh PetscBool ismpikaij,isseqkaij; 5149bd79ccSDebojyoti Ghosh 5249bd79ccSDebojyoti Ghosh PetscFunctionBegin; 5349bd79ccSDebojyoti Ghosh ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij);CHKERRQ(ierr); 5449bd79ccSDebojyoti Ghosh ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQKAIJ,&isseqkaij);CHKERRQ(ierr); 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"); 6449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 6549bd79ccSDebojyoti Ghosh } 6649bd79ccSDebojyoti Ghosh 6749bd79ccSDebojyoti Ghosh /*@C 6849bd79ccSDebojyoti Ghosh MatKAIJGetS - Get the S matrix describing the shift action of the KAIJ 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: 7349bd79ccSDebojyoti Ghosh . A - the KAIJ 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 80a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 8131a97b9aSRichard Tran Mills 8249bd79ccSDebojyoti Ghosh Level: advanced 8349bd79ccSDebojyoti Ghosh 8431a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 8549bd79ccSDebojyoti Ghosh @*/ 86a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetS(Mat A,PetscInt *m,PetscInt *n,PetscScalar **S) 8749bd79ccSDebojyoti Ghosh { 8849bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 8949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 90a5b5c723SRichard Tran Mills if (m) *m = b->p; 91a5b5c723SRichard Tran Mills if (n) *n = b->q; 92a5b5c723SRichard Tran Mills if (S) *S = b->S; 93a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 94a5b5c723SRichard Tran Mills } 95a5b5c723SRichard Tran Mills 96a5b5c723SRichard Tran Mills /*@C 97a5b5c723SRichard Tran Mills MatKAIJGetSRead - Get a read-only pointer to the S matrix describing the shift action of the KAIJ matrix 98a5b5c723SRichard Tran Mills 99a5b5c723SRichard Tran Mills Not Collective; the entire S is stored and returned independently on all processes. 100a5b5c723SRichard Tran Mills 101a5b5c723SRichard Tran Mills Input Parameter: 102a5b5c723SRichard Tran Mills . A - the KAIJ matrix 103a5b5c723SRichard Tran Mills 104a5b5c723SRichard Tran Mills Output Parameters: 105a5b5c723SRichard Tran Mills + m - the number of rows in S 106a5b5c723SRichard Tran Mills . n - the number of columns in S 107a5b5c723SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 108a5b5c723SRichard Tran Mills 109a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 110a5b5c723SRichard Tran Mills 111a5b5c723SRichard Tran Mills Level: advanced 112a5b5c723SRichard Tran Mills 113a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 114a5b5c723SRichard Tran Mills @*/ 115a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetSRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **S) 116a5b5c723SRichard Tran Mills { 117a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 118a5b5c723SRichard Tran Mills PetscFunctionBegin; 119a5b5c723SRichard Tran Mills if (m) *m = b->p; 120a5b5c723SRichard Tran Mills if (n) *n = b->q; 121a5b5c723SRichard Tran Mills if (S) *S = b->S; 122a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 123a5b5c723SRichard Tran Mills } 124a5b5c723SRichard Tran Mills 125a5b5c723SRichard Tran Mills /*@C 126a5b5c723SRichard Tran Mills MatKAIJRestoreS - Restore array obtained with MatKAIJGetS() 127a5b5c723SRichard Tran Mills 128a5b5c723SRichard Tran Mills Not collective 129a5b5c723SRichard Tran Mills 130a5b5c723SRichard Tran Mills Input Parameter: 131a5b5c723SRichard Tran Mills . A - the KAIJ matrix 132a5b5c723SRichard Tran Mills 133a5b5c723SRichard Tran Mills Output Parameter: 134a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS() 135a5b5c723SRichard Tran Mills 136a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 137a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 138a5b5c723SRichard Tran Mills 139a5b5c723SRichard Tran Mills Level: advanced 140a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead() 141a5b5c723SRichard Tran Mills @*/ 142a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreS(Mat A,PetscScalar **S) 143a5b5c723SRichard Tran Mills { 14466f58c76SRichard Tran Mills PetscErrorCode ierr; 14566f58c76SRichard Tran Mills 146a5b5c723SRichard Tran Mills PetscFunctionBegin; 147a5b5c723SRichard Tran Mills if (S) *S = NULL; 14866f58c76SRichard Tran Mills ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 149a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 150a5b5c723SRichard Tran Mills } 151a5b5c723SRichard Tran Mills 152a5b5c723SRichard Tran Mills /*@C 153a5b5c723SRichard Tran Mills MatKAIJRestoreSRead - Restore array obtained with MatKAIJGetSRead() 154a5b5c723SRichard Tran Mills 155a5b5c723SRichard Tran Mills Not collective 156a5b5c723SRichard Tran Mills 157a5b5c723SRichard Tran Mills Input Parameter: 158a5b5c723SRichard Tran Mills . A - the KAIJ matrix 159a5b5c723SRichard Tran Mills 160a5b5c723SRichard Tran Mills Output Parameter: 161a5b5c723SRichard Tran Mills . S - location of pointer to array obtained with MatKAIJGetS() 162a5b5c723SRichard Tran Mills 163a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 164a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 165a5b5c723SRichard Tran Mills 166a5b5c723SRichard Tran Mills Level: advanced 167a5b5c723SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJGetSRead(), MatKAIJRestoreSRead() 168a5b5c723SRichard Tran Mills @*/ 169a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreSRead(Mat A,const PetscScalar **S) 170a5b5c723SRichard Tran Mills { 171a5b5c723SRichard Tran Mills PetscFunctionBegin; 172a5b5c723SRichard Tran Mills if (S) *S = NULL; 17349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 17449bd79ccSDebojyoti Ghosh } 17549bd79ccSDebojyoti Ghosh 17649bd79ccSDebojyoti Ghosh /*@C 17731a97b9aSRichard Tran Mills MatKAIJGetT - Get the transformation matrix T associated with the KAIJ matrix 17849bd79ccSDebojyoti Ghosh 1790567c835SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 18049bd79ccSDebojyoti Ghosh 18149bd79ccSDebojyoti Ghosh Input Parameter: 18249bd79ccSDebojyoti Ghosh . A - the KAIJ matrix 18349bd79ccSDebojyoti Ghosh 184d8d19677SJose E. Roman Output Parameters: 185a5b5c723SRichard Tran Mills + m - the number of rows in T 186a5b5c723SRichard Tran Mills . n - the number of columns in T 187a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 18849bd79ccSDebojyoti Ghosh 189a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 19031a97b9aSRichard Tran Mills 19149bd79ccSDebojyoti Ghosh Level: advanced 19249bd79ccSDebojyoti Ghosh 19331a97b9aSRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 19449bd79ccSDebojyoti Ghosh @*/ 195a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetT(Mat A,PetscInt *m,PetscInt *n,PetscScalar **T) 19649bd79ccSDebojyoti Ghosh { 19749bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 19849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 199a5b5c723SRichard Tran Mills if (m) *m = b->p; 200a5b5c723SRichard Tran Mills if (n) *n = b->q; 201a5b5c723SRichard Tran Mills if (T) *T = b->T; 202a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 203a5b5c723SRichard Tran Mills } 204a5b5c723SRichard Tran Mills 205a5b5c723SRichard Tran Mills /*@C 206a5b5c723SRichard Tran Mills MatKAIJGetTRead - Get a read-only pointer to the transformation matrix T associated with the KAIJ matrix 207a5b5c723SRichard Tran Mills 208a5b5c723SRichard Tran Mills Not Collective; the entire T is stored and returned independently on all processes 209a5b5c723SRichard Tran Mills 210a5b5c723SRichard Tran Mills Input Parameter: 211a5b5c723SRichard Tran Mills . A - the KAIJ matrix 212a5b5c723SRichard Tran Mills 213d8d19677SJose E. Roman Output Parameters: 214a5b5c723SRichard Tran Mills + m - the number of rows in T 215a5b5c723SRichard Tran Mills . n - the number of columns in T 216a5b5c723SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 217a5b5c723SRichard Tran Mills 218a5b5c723SRichard Tran Mills Note: All output parameters are optional (pass NULL or PETSC_IGNORE if not desired) 219a5b5c723SRichard Tran Mills 220a5b5c723SRichard Tran Mills Level: advanced 221a5b5c723SRichard Tran Mills 222a5b5c723SRichard Tran Mills .seealso: MatCreateKAIJ(), MatGetBlockSizes() 223a5b5c723SRichard Tran Mills @*/ 224a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJGetTRead(Mat A,PetscInt *m,PetscInt *n,const PetscScalar **T) 225a5b5c723SRichard Tran Mills { 226a5b5c723SRichard Tran Mills Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 227a5b5c723SRichard Tran Mills PetscFunctionBegin; 228a5b5c723SRichard Tran Mills if (m) *m = b->p; 229a5b5c723SRichard Tran Mills if (n) *n = b->q; 230a5b5c723SRichard Tran Mills if (T) *T = b->T; 231a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 232a5b5c723SRichard Tran Mills } 233a5b5c723SRichard Tran Mills 234a5b5c723SRichard Tran Mills /*@C 235a5b5c723SRichard Tran Mills MatKAIJRestoreT - Restore array obtained with MatKAIJGetT() 236a5b5c723SRichard Tran Mills 237a5b5c723SRichard Tran Mills Not collective 238a5b5c723SRichard Tran Mills 239a5b5c723SRichard Tran Mills Input Parameter: 240a5b5c723SRichard Tran Mills . A - the KAIJ matrix 241a5b5c723SRichard Tran Mills 242a5b5c723SRichard Tran Mills Output Parameter: 243a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS() 244a5b5c723SRichard Tran Mills 245a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 246a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 247a5b5c723SRichard Tran Mills 248a5b5c723SRichard Tran Mills Level: advanced 249a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead() 250a5b5c723SRichard Tran Mills @*/ 251a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreT(Mat A,PetscScalar **T) 252a5b5c723SRichard Tran Mills { 25366f58c76SRichard Tran Mills PetscErrorCode ierr; 25466f58c76SRichard Tran Mills 255a5b5c723SRichard Tran Mills PetscFunctionBegin; 256a5b5c723SRichard Tran Mills if (T) *T = NULL; 25766f58c76SRichard Tran Mills ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 258a5b5c723SRichard Tran Mills PetscFunctionReturn(0); 259a5b5c723SRichard Tran Mills } 260a5b5c723SRichard Tran Mills 261a5b5c723SRichard Tran Mills /*@C 262a5b5c723SRichard Tran Mills MatKAIJRestoreTRead - Restore array obtained with MatKAIJGetTRead() 263a5b5c723SRichard Tran Mills 264a5b5c723SRichard Tran Mills Not collective 265a5b5c723SRichard Tran Mills 266a5b5c723SRichard Tran Mills Input Parameter: 267a5b5c723SRichard Tran Mills . A - the KAIJ matrix 268a5b5c723SRichard Tran Mills 269a5b5c723SRichard Tran Mills Output Parameter: 270a5b5c723SRichard Tran Mills . T - location of pointer to array obtained with MatKAIJGetS() 271a5b5c723SRichard Tran Mills 272a5b5c723SRichard Tran Mills Note: This routine zeros the array pointer to prevent accidental reuse after it has been restored. 273a5b5c723SRichard Tran Mills If NULL is passed, it will not attempt to zero the array pointer. 274a5b5c723SRichard Tran Mills 275a5b5c723SRichard Tran Mills Level: advanced 276a5b5c723SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJGetTRead(), MatKAIJRestoreTRead() 277a5b5c723SRichard Tran Mills @*/ 278a5b5c723SRichard Tran Mills PetscErrorCode MatKAIJRestoreTRead(Mat A,const PetscScalar **T) 279a5b5c723SRichard Tran Mills { 280a5b5c723SRichard Tran Mills PetscFunctionBegin; 281a5b5c723SRichard Tran Mills if (T) *T = NULL; 28249bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 28349bd79ccSDebojyoti Ghosh } 28449bd79ccSDebojyoti Ghosh 2850567c835SRichard Tran Mills /*@ 2860567c835SRichard Tran Mills MatKAIJSetAIJ - Set the AIJ matrix describing the blockwise action of the KAIJ matrix 2870567c835SRichard Tran Mills 2880567c835SRichard Tran Mills Logically Collective; if the AIJ matrix is parallel, the KAIJ matrix is also parallel 2890567c835SRichard Tran Mills 2900567c835SRichard Tran Mills Input Parameters: 2910567c835SRichard Tran Mills + A - the KAIJ matrix 2920567c835SRichard Tran Mills - B - the AIJ matrix 2930567c835SRichard Tran Mills 29415b9d025SRichard Tran Mills Notes: 29515b9d025SRichard 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. 29615b9d025SRichard Tran Mills Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix. 29715b9d025SRichard Tran Mills 2980567c835SRichard Tran Mills Level: advanced 2990567c835SRichard Tran Mills 3000567c835SRichard Tran Mills .seealso: MatKAIJGetAIJ(), MatKAIJSetS(), MatKAIJSetT() 3010567c835SRichard Tran Mills @*/ 3020567c835SRichard Tran Mills PetscErrorCode MatKAIJSetAIJ(Mat A,Mat B) 3030567c835SRichard Tran Mills { 3040567c835SRichard Tran Mills PetscErrorCode ierr; 3050567c835SRichard Tran Mills PetscMPIInt size; 306*06d61a37SPierre Jolivet PetscBool flg; 3070567c835SRichard Tran Mills 3080567c835SRichard Tran Mills PetscFunctionBegin; 309ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 3100567c835SRichard Tran Mills if (size == 1) { 311*06d61a37SPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&flg);CHKERRQ(ierr); 312*06d61a37SPierre Jolivet if (!flg) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MatKAIJSetAIJ() with MATSEQKAIJ does not support %s as the AIJ mat",((PetscObject)B)->type_name); 3130567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 3140567c835SRichard Tran Mills a->AIJ = B; 3150567c835SRichard Tran Mills } else { 3160567c835SRichard Tran Mills Mat_MPIKAIJ *a = (Mat_MPIKAIJ*)A->data; 3170567c835SRichard Tran Mills a->A = B; 3180567c835SRichard Tran Mills } 31915b9d025SRichard Tran Mills ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 3200567c835SRichard Tran Mills PetscFunctionReturn(0); 3210567c835SRichard Tran Mills } 3220567c835SRichard Tran Mills 3230567c835SRichard Tran Mills /*@C 3240567c835SRichard Tran Mills MatKAIJSetS - Set the S matrix describing the shift action of the KAIJ matrix 3250567c835SRichard Tran Mills 3260567c835SRichard Tran Mills Logically Collective; the entire S is stored independently on all processes. 3270567c835SRichard Tran Mills 3280567c835SRichard Tran Mills Input Parameters: 3290567c835SRichard Tran Mills + A - the KAIJ matrix 3300567c835SRichard Tran Mills . p - the number of rows in S 3310567c835SRichard Tran Mills . q - the number of columns in S 3320567c835SRichard Tran Mills - S - the S matrix, in form of a scalar array in column-major format 3330567c835SRichard Tran Mills 3340567c835SRichard Tran Mills Notes: The dimensions p and q must match those of the transformation matrix T associated with the KAIJ matrix. 33588f48298SRichard Tran Mills The S matrix is copied, so the user can destroy this array. 3360567c835SRichard Tran Mills 3370567c835SRichard Tran Mills Level: Advanced 3380567c835SRichard Tran Mills 3390567c835SRichard Tran Mills .seealso: MatKAIJGetS(), MatKAIJSetT(), MatKAIJSetAIJ() 3400567c835SRichard Tran Mills @*/ 3410567c835SRichard Tran Mills PetscErrorCode MatKAIJSetS(Mat A,PetscInt p,PetscInt q,const PetscScalar S[]) 3420567c835SRichard Tran Mills { 3430567c835SRichard Tran Mills PetscErrorCode ierr; 3440567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 3450567c835SRichard Tran Mills 3460567c835SRichard Tran Mills PetscFunctionBegin; 3470567c835SRichard Tran Mills ierr = PetscFree(a->S);CHKERRQ(ierr); 3480567c835SRichard Tran Mills if (S) { 349a84f8069SRichard Tran Mills ierr = PetscMalloc1(p*q*sizeof(PetscScalar),&a->S);CHKERRQ(ierr); 3500567c835SRichard Tran Mills ierr = PetscMemcpy(a->S,S,p*q*sizeof(PetscScalar));CHKERRQ(ierr); 3510567c835SRichard Tran Mills } else a->S = NULL; 3520567c835SRichard Tran Mills 3530567c835SRichard Tran Mills a->p = p; 3540567c835SRichard Tran Mills a->q = q; 3550567c835SRichard Tran Mills PetscFunctionReturn(0); 3560567c835SRichard Tran Mills } 3570567c835SRichard Tran Mills 3580567c835SRichard Tran Mills /*@C 359910cf402Sprj- MatKAIJGetScaledIdentity - Check if both S and T are scaled identities. 360910cf402Sprj- 361910cf402Sprj- Logically Collective. 362910cf402Sprj- 363910cf402Sprj- Input Parameter: 364910cf402Sprj- . A - the KAIJ matrix 365910cf402Sprj- 366910cf402Sprj- Output Parameter: 367910cf402Sprj- . identity - the Boolean value 368910cf402Sprj- 369910cf402Sprj- Level: Advanced 370910cf402Sprj- 371910cf402Sprj- .seealso: MatKAIJGetS(), MatKAIJGetT() 372910cf402Sprj- @*/ 373910cf402Sprj- PetscErrorCode MatKAIJGetScaledIdentity(Mat A,PetscBool* identity) 374910cf402Sprj- { 375910cf402Sprj- Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 376910cf402Sprj- PetscInt i,j; 377910cf402Sprj- 378910cf402Sprj- PetscFunctionBegin; 379910cf402Sprj- if (a->p != a->q) { 380910cf402Sprj- *identity = PETSC_FALSE; 381910cf402Sprj- PetscFunctionReturn(0); 382910cf402Sprj- } else *identity = PETSC_TRUE; 383910cf402Sprj- if (!a->isTI || a->S) { 384910cf402Sprj- for (i=0; i<a->p && *identity; i++) { 385910cf402Sprj- for (j=0; j<a->p && *identity; j++) { 386910cf402Sprj- if (i != j) { 387910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 388910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i+j*a->p]) > PETSC_SMALL) *identity = PETSC_FALSE; 389910cf402Sprj- } else { 390910cf402Sprj- if (a->S && PetscAbsScalar(a->S[i*(a->p+1)]-a->S[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 391910cf402Sprj- if (a->T && PetscAbsScalar(a->T[i*(a->p+1)]-a->T[0]) > PETSC_SMALL) *identity = PETSC_FALSE; 392910cf402Sprj- } 393910cf402Sprj- } 394910cf402Sprj- } 395910cf402Sprj- } 396910cf402Sprj- PetscFunctionReturn(0); 397910cf402Sprj- } 398910cf402Sprj- 399910cf402Sprj- /*@C 4000567c835SRichard Tran Mills MatKAIJSetT - Set the transformation matrix T associated with the KAIJ matrix 4010567c835SRichard Tran Mills 4020567c835SRichard Tran Mills Logically Collective; the entire T is stored independently on all processes. 4030567c835SRichard Tran Mills 4040567c835SRichard Tran Mills Input Parameters: 4050567c835SRichard Tran Mills + A - the KAIJ matrix 4060567c835SRichard Tran Mills . p - the number of rows in S 4070567c835SRichard Tran Mills . q - the number of columns in S 4080567c835SRichard Tran Mills - T - the T matrix, in form of a scalar array in column-major format 4090567c835SRichard Tran Mills 4100567c835SRichard Tran Mills Notes: The dimensions p and q must match those of the shift matrix S associated with the KAIJ matrix. 41188f48298SRichard Tran Mills The T matrix is copied, so the user can destroy this array. 4120567c835SRichard Tran Mills 4130567c835SRichard Tran Mills Level: Advanced 4140567c835SRichard Tran Mills 4150567c835SRichard Tran Mills .seealso: MatKAIJGetT(), MatKAIJSetS(), MatKAIJSetAIJ() 4160567c835SRichard Tran Mills @*/ 4170567c835SRichard Tran Mills PetscErrorCode MatKAIJSetT(Mat A,PetscInt p,PetscInt q,const PetscScalar T[]) 4180567c835SRichard Tran Mills { 4190567c835SRichard Tran Mills PetscErrorCode ierr; 4200567c835SRichard Tran Mills PetscInt i,j; 4210567c835SRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 4220567c835SRichard Tran Mills PetscBool isTI = PETSC_FALSE; 4230567c835SRichard Tran Mills 4240567c835SRichard Tran Mills PetscFunctionBegin; 4250567c835SRichard Tran Mills /* check if T is an identity matrix */ 4260567c835SRichard Tran Mills if (T && (p == q)) { 4270567c835SRichard Tran Mills isTI = PETSC_TRUE; 4280567c835SRichard Tran Mills for (i=0; i<p; i++) { 4290567c835SRichard Tran Mills for (j=0; j<q; j++) { 4300567c835SRichard Tran Mills if (i == j) { 4310567c835SRichard Tran Mills /* diagonal term must be 1 */ 4320567c835SRichard Tran Mills if (T[i+j*p] != 1.0) isTI = PETSC_FALSE; 4330567c835SRichard Tran Mills } else { 4340567c835SRichard Tran Mills /* off-diagonal term must be 0 */ 4350567c835SRichard Tran Mills if (T[i+j*p] != 0.0) isTI = PETSC_FALSE; 4360567c835SRichard Tran Mills } 4370567c835SRichard Tran Mills } 4380567c835SRichard Tran Mills } 4390567c835SRichard Tran Mills } 4400567c835SRichard Tran Mills a->isTI = isTI; 4410567c835SRichard Tran Mills 4420567c835SRichard Tran Mills ierr = PetscFree(a->T);CHKERRQ(ierr); 4430567c835SRichard Tran Mills if (T && (!isTI)) { 444a84f8069SRichard Tran Mills ierr = PetscMalloc1(p*q*sizeof(PetscScalar),&a->T);CHKERRQ(ierr); 4450567c835SRichard Tran Mills ierr = PetscMemcpy(a->T,T,p*q*sizeof(PetscScalar));CHKERRQ(ierr); 44650d19d74SRichard Tran Mills } else a->T = NULL; 4470567c835SRichard Tran Mills 4480567c835SRichard Tran Mills a->p = p; 4490567c835SRichard Tran Mills a->q = q; 4500567c835SRichard Tran Mills PetscFunctionReturn(0); 4510567c835SRichard Tran Mills } 4520567c835SRichard Tran Mills 45349bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_SeqKAIJ(Mat A) 45449bd79ccSDebojyoti Ghosh { 45549bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 45649bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 45749bd79ccSDebojyoti Ghosh 45849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 45949bd79ccSDebojyoti Ghosh ierr = MatDestroy(&b->AIJ);CHKERRQ(ierr); 460a84f8069SRichard Tran Mills ierr = PetscFree(b->S);CHKERRQ(ierr); 461a84f8069SRichard Tran Mills ierr = PetscFree(b->T);CHKERRQ(ierr); 462a84f8069SRichard Tran Mills ierr = PetscFree(b->ibdiag);CHKERRQ(ierr); 46349bd79ccSDebojyoti Ghosh ierr = PetscFree5(b->sor.w,b->sor.y,b->sor.work,b->sor.t,b->sor.arr);CHKERRQ(ierr); 464*06d61a37SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqkaij_seqaij_C",NULL);CHKERRQ(ierr); 46549bd79ccSDebojyoti Ghosh ierr = PetscFree(A->data);CHKERRQ(ierr); 46649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 46749bd79ccSDebojyoti Ghosh } 46849bd79ccSDebojyoti Ghosh 469e0e5a793SRichard Tran Mills PETSC_INTERN PetscErrorCode MatKAIJ_build_AIJ_OAIJ(Mat A) 470e0e5a793SRichard Tran Mills { 471e0e5a793SRichard Tran Mills PetscErrorCode ierr; 472e0e5a793SRichard Tran Mills Mat_MPIKAIJ *a; 473e0e5a793SRichard Tran Mills Mat_MPIAIJ *mpiaij; 474e0e5a793SRichard Tran Mills PetscScalar *T; 475e0e5a793SRichard Tran Mills PetscInt i,j; 476e0e5a793SRichard Tran Mills PetscObjectState state; 477e0e5a793SRichard Tran Mills 478e0e5a793SRichard Tran Mills PetscFunctionBegin; 479e0e5a793SRichard Tran Mills a = (Mat_MPIKAIJ*)A->data; 480e0e5a793SRichard Tran Mills mpiaij = (Mat_MPIAIJ*)a->A->data; 481e0e5a793SRichard Tran Mills 482e0e5a793SRichard Tran Mills ierr = PetscObjectStateGet((PetscObject)a->A,&state);CHKERRQ(ierr); 483e0e5a793SRichard Tran Mills if (state == a->state) { 484e0e5a793SRichard Tran Mills /* The existing AIJ and KAIJ members are up-to-date, so simply exit. */ 485e0e5a793SRichard Tran Mills PetscFunctionReturn(0); 486e0e5a793SRichard Tran Mills } else { 487e0e5a793SRichard Tran Mills ierr = MatDestroy(&a->AIJ);CHKERRQ(ierr); 488e0e5a793SRichard Tran Mills ierr = MatDestroy(&a->OAIJ);CHKERRQ(ierr); 489e0e5a793SRichard Tran Mills if (a->isTI) { 490e0e5a793SRichard Tran Mills /* If the transformation matrix associated with the parallel matrix A is the identity matrix, then a->T will be NULL. 491e0e5a793SRichard Tran Mills * In this case, if we pass a->T directly to the MatCreateKAIJ() calls to create the sequential submatrices, the routine will 492e0e5a793SRichard Tran Mills * not be able to tell that transformation matrix should be set to the identity; thus we create a temporary identity matrix 493e0e5a793SRichard Tran Mills * to pass in. */ 494e0e5a793SRichard Tran Mills ierr = PetscMalloc1(a->p*a->q*sizeof(PetscScalar),&T);CHKERRQ(ierr); 495e0e5a793SRichard Tran Mills for (i=0; i<a->p; i++) { 496e0e5a793SRichard Tran Mills for (j=0; j<a->q; j++) { 497e0e5a793SRichard Tran Mills if (i==j) T[i+j*a->p] = 1.0; 498e0e5a793SRichard Tran Mills else T[i+j*a->p] = 0.0; 499e0e5a793SRichard Tran Mills } 500e0e5a793SRichard Tran Mills } 501e0e5a793SRichard Tran Mills } else T = a->T; 502e0e5a793SRichard Tran Mills ierr = MatCreateKAIJ(mpiaij->A,a->p,a->q,a->S,T,&a->AIJ);CHKERRQ(ierr); 503e0e5a793SRichard Tran Mills ierr = MatCreateKAIJ(mpiaij->B,a->p,a->q,NULL,T,&a->OAIJ);CHKERRQ(ierr); 504e0e5a793SRichard Tran Mills if (a->isTI) { 505e0e5a793SRichard Tran Mills ierr = PetscFree(T);CHKERRQ(ierr); 506e0e5a793SRichard Tran Mills } 507e0e5a793SRichard Tran Mills a->state = state; 508e0e5a793SRichard Tran Mills } 509e0e5a793SRichard Tran Mills 510e0e5a793SRichard Tran Mills PetscFunctionReturn(0); 511e0e5a793SRichard Tran Mills } 512e0e5a793SRichard Tran Mills 51349bd79ccSDebojyoti Ghosh PetscErrorCode MatSetUp_KAIJ(Mat A) 51449bd79ccSDebojyoti Ghosh { 5150567c835SRichard Tran Mills PetscErrorCode ierr; 5160567c835SRichard Tran Mills PetscInt n; 5170567c835SRichard Tran Mills PetscMPIInt size; 5180567c835SRichard Tran Mills Mat_SeqKAIJ *seqkaij = (Mat_SeqKAIJ*)A->data; 5190567c835SRichard Tran Mills 52049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 521ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 5220567c835SRichard Tran Mills if (size == 1) { 5230567c835SRichard Tran Mills ierr = 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);CHKERRQ(ierr); 5240567c835SRichard Tran Mills ierr = PetscLayoutSetBlockSize(A->rmap,seqkaij->p);CHKERRQ(ierr); 5250567c835SRichard Tran Mills ierr = PetscLayoutSetBlockSize(A->cmap,seqkaij->q);CHKERRQ(ierr); 5260567c835SRichard Tran Mills ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5270567c835SRichard Tran Mills ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5280567c835SRichard Tran Mills } else { 5290567c835SRichard Tran Mills Mat_MPIKAIJ *a; 5300567c835SRichard Tran Mills Mat_MPIAIJ *mpiaij; 5310567c835SRichard Tran Mills IS from,to; 5320567c835SRichard Tran Mills Vec gvec; 5330567c835SRichard Tran Mills 5340567c835SRichard Tran Mills a = (Mat_MPIKAIJ*)A->data; 535d3f912faSRichard Tran Mills mpiaij = (Mat_MPIAIJ*)a->A->data; 5360567c835SRichard Tran Mills ierr = 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);CHKERRQ(ierr); 5370567c835SRichard Tran Mills ierr = PetscLayoutSetBlockSize(A->rmap,seqkaij->p);CHKERRQ(ierr); 5380567c835SRichard Tran Mills ierr = PetscLayoutSetBlockSize(A->cmap,seqkaij->q);CHKERRQ(ierr); 5390567c835SRichard Tran Mills ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5400567c835SRichard Tran Mills ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5410567c835SRichard Tran Mills 542eacc1826SRichard Tran Mills ierr = MatKAIJ_build_AIJ_OAIJ(A);CHKERRQ(ierr); 5430567c835SRichard Tran Mills 5440567c835SRichard Tran Mills ierr = VecGetSize(mpiaij->lvec,&n);CHKERRQ(ierr); 5450567c835SRichard Tran Mills ierr = VecCreate(PETSC_COMM_SELF,&a->w);CHKERRQ(ierr); 5460567c835SRichard Tran Mills ierr = VecSetSizes(a->w,n*a->q,n*a->q);CHKERRQ(ierr); 5470567c835SRichard Tran Mills ierr = VecSetBlockSize(a->w,a->q);CHKERRQ(ierr); 5480567c835SRichard Tran Mills ierr = VecSetType(a->w,VECSEQ);CHKERRQ(ierr); 5490567c835SRichard Tran Mills 5500567c835SRichard Tran Mills /* create two temporary Index sets for build scatter gather */ 5510567c835SRichard Tran Mills ierr = ISCreateBlock(PetscObjectComm((PetscObject)a->A),a->q,n,mpiaij->garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 5520567c835SRichard Tran Mills ierr = ISCreateStride(PETSC_COMM_SELF,n*a->q,0,1,&to);CHKERRQ(ierr); 5530567c835SRichard Tran Mills 5540567c835SRichard Tran Mills /* create temporary global vector to generate scatter context */ 5550567c835SRichard Tran Mills ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)a->A),a->q,a->q*a->A->cmap->n,a->q*a->A->cmap->N,NULL,&gvec);CHKERRQ(ierr); 5560567c835SRichard Tran Mills 5570567c835SRichard Tran Mills /* generate the scatter context */ 5584589b4e5SRichard Tran Mills ierr = VecScatterCreate(gvec,from,a->w,to,&a->ctx);CHKERRQ(ierr); 5590567c835SRichard Tran Mills 5600567c835SRichard Tran Mills ierr = ISDestroy(&from);CHKERRQ(ierr); 5610567c835SRichard Tran Mills ierr = ISDestroy(&to);CHKERRQ(ierr); 5620567c835SRichard Tran Mills ierr = VecDestroy(&gvec);CHKERRQ(ierr); 5630567c835SRichard Tran Mills } 5640567c835SRichard Tran Mills 5650567c835SRichard Tran Mills A->assembled = PETSC_TRUE; 56649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 56749bd79ccSDebojyoti Ghosh } 56849bd79ccSDebojyoti Ghosh 569e6985dafSRichard Tran Mills PetscErrorCode MatView_KAIJ(Mat A,PetscViewer viewer) 57049bd79ccSDebojyoti Ghosh { 571e6985dafSRichard Tran Mills PetscViewerFormat format; 572e6985dafSRichard Tran Mills Mat_SeqKAIJ *a = (Mat_SeqKAIJ*)A->data; 57349bd79ccSDebojyoti Ghosh Mat B; 574e6985dafSRichard Tran Mills PetscInt i; 575e6985dafSRichard Tran Mills PetscErrorCode ierr; 576e6985dafSRichard Tran Mills PetscBool ismpikaij; 57749bd79ccSDebojyoti Ghosh 57849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 579e6985dafSRichard Tran Mills ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIKAIJ,&ismpikaij);CHKERRQ(ierr); 580e6985dafSRichard Tran Mills ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 581e6985dafSRichard Tran Mills if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_IMPL) { 582c0aa6a63SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(viewer,"S and T have %" PetscInt_FMT " rows and %" PetscInt_FMT " columns\n",a->p,a->q);CHKERRQ(ierr); 583e6985dafSRichard Tran Mills 584e6985dafSRichard Tran Mills /* Print appropriate details for S. */ 585e6985dafSRichard Tran Mills if (!a->S) { 5862ae760e3SRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"S is NULL\n");CHKERRQ(ierr); 587e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 588e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"Entries of S are ");CHKERRQ(ierr); 589e6985dafSRichard Tran Mills for (i=0; i<(a->p * a->q); i++) { 590e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 591e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->S[i]),(double)PetscImaginaryPart(a->S[i]));CHKERRQ(ierr); 592e6985dafSRichard Tran Mills #else 593e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->S[i]));CHKERRQ(ierr); 594e6985dafSRichard Tran Mills #endif 595e6985dafSRichard Tran Mills } 596e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 59749bd79ccSDebojyoti Ghosh } 59849bd79ccSDebojyoti Ghosh 599e6985dafSRichard Tran Mills /* Print appropriate details for T. */ 600e6985dafSRichard Tran Mills if (a->isTI) { 6012ae760e3SRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"T is the identity matrix\n");CHKERRQ(ierr); 602e6985dafSRichard Tran Mills } else if (!a->T) { 6032ae760e3SRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"T is NULL\n");CHKERRQ(ierr); 604e6985dafSRichard Tran Mills } else if (format == PETSC_VIEWER_ASCII_IMPL) { 605e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"Entries of T are ");CHKERRQ(ierr); 606e6985dafSRichard Tran Mills for (i=0; i<(a->p * a->q); i++) { 607e6985dafSRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 608e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"%18.16e %18.16e ",(double)PetscRealPart(a->T[i]),(double)PetscImaginaryPart(a->T[i]));CHKERRQ(ierr); 609e6985dafSRichard Tran Mills #else 610e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(a->T[i]));CHKERRQ(ierr); 611e6985dafSRichard Tran Mills #endif 612e6985dafSRichard Tran Mills } 613e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 614e6985dafSRichard Tran Mills } 61549bd79ccSDebojyoti Ghosh 616e6985dafSRichard Tran Mills /* Now print details for the AIJ matrix, using the AIJ viewer. */ 617e6985dafSRichard Tran Mills ierr = PetscViewerASCIIPrintf(viewer,"Now viewing the associated AIJ matrix:\n");CHKERRQ(ierr); 618e6985dafSRichard Tran Mills if (ismpikaij) { 619e6985dafSRichard Tran Mills Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 620e6985dafSRichard Tran Mills ierr = MatView(b->A,viewer);CHKERRQ(ierr); 621e6985dafSRichard Tran Mills } else { 622e6985dafSRichard Tran Mills ierr = MatView(a->AIJ,viewer);CHKERRQ(ierr); 623e6985dafSRichard Tran Mills } 624e6985dafSRichard Tran Mills 625e6985dafSRichard Tran Mills } else { 626e6985dafSRichard Tran Mills /* For all other matrix viewer output formats, simply convert to an AIJ matrix and call MatView() on that. */ 627*06d61a37SPierre Jolivet ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 62849bd79ccSDebojyoti Ghosh ierr = MatView(B,viewer);CHKERRQ(ierr); 62949bd79ccSDebojyoti Ghosh ierr = MatDestroy(&B);CHKERRQ(ierr); 630e6985dafSRichard Tran Mills } 63149bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 63249bd79ccSDebojyoti Ghosh } 63349bd79ccSDebojyoti Ghosh 63449bd79ccSDebojyoti Ghosh PetscErrorCode MatDestroy_MPIKAIJ(Mat A) 63549bd79ccSDebojyoti Ghosh { 63649bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 63749bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 63849bd79ccSDebojyoti Ghosh 63949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 64049bd79ccSDebojyoti Ghosh ierr = MatDestroy(&b->AIJ);CHKERRQ(ierr); 64149bd79ccSDebojyoti Ghosh ierr = MatDestroy(&b->OAIJ);CHKERRQ(ierr); 64249bd79ccSDebojyoti Ghosh ierr = MatDestroy(&b->A);CHKERRQ(ierr); 64349bd79ccSDebojyoti Ghosh ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 64449bd79ccSDebojyoti Ghosh ierr = VecDestroy(&b->w);CHKERRQ(ierr); 645a84f8069SRichard Tran Mills ierr = PetscFree(b->S);CHKERRQ(ierr); 646a84f8069SRichard Tran Mills ierr = PetscFree(b->T);CHKERRQ(ierr); 647a84f8069SRichard Tran Mills ierr = PetscFree(b->ibdiag);CHKERRQ(ierr); 648*06d61a37SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetDiagonalBlock_C",NULL);CHKERRQ(ierr); 649*06d61a37SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_mpikaij_mpiaij_C",NULL);CHKERRQ(ierr); 65049bd79ccSDebojyoti Ghosh ierr = PetscFree(A->data);CHKERRQ(ierr); 65149bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 65249bd79ccSDebojyoti Ghosh } 65349bd79ccSDebojyoti Ghosh 65449bd79ccSDebojyoti Ghosh /* --------------------------------------------------------------------------------------*/ 65549bd79ccSDebojyoti Ghosh 65649bd79ccSDebojyoti Ghosh /* zz = yy + Axx */ 657836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_SeqKAIJ(Mat A,Vec xx,Vec yy,Vec zz) 65849bd79ccSDebojyoti Ghosh { 65949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 66049bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 66149bd79ccSDebojyoti Ghosh const PetscScalar *s = b->S, *t = b->T; 66249bd79ccSDebojyoti Ghosh const PetscScalar *x,*v,*bx; 66349bd79ccSDebojyoti Ghosh PetscScalar *y,*sums; 66449bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 66549bd79ccSDebojyoti Ghosh const PetscInt m = b->AIJ->rmap->n,*idx,*ii; 66649bd79ccSDebojyoti Ghosh PetscInt n,i,jrow,j,l,p=b->p,q=b->q,k; 66749bd79ccSDebojyoti Ghosh 66849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 66949bd79ccSDebojyoti Ghosh if (!yy) { 67049bd79ccSDebojyoti Ghosh ierr = VecSet(zz,0.0);CHKERRQ(ierr); 67149bd79ccSDebojyoti Ghosh } else { 67249bd79ccSDebojyoti Ghosh ierr = VecCopy(yy,zz);CHKERRQ(ierr); 67349bd79ccSDebojyoti Ghosh } 67449bd79ccSDebojyoti Ghosh if ((!s) && (!t) && (!b->isTI)) PetscFunctionReturn(0); 67549bd79ccSDebojyoti Ghosh 67649bd79ccSDebojyoti Ghosh ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 67749bd79ccSDebojyoti Ghosh ierr = VecGetArray(zz,&y);CHKERRQ(ierr); 67849bd79ccSDebojyoti Ghosh idx = a->j; 67949bd79ccSDebojyoti Ghosh v = a->a; 68049bd79ccSDebojyoti Ghosh ii = a->i; 68149bd79ccSDebojyoti Ghosh 68249bd79ccSDebojyoti Ghosh if (b->isTI) { 68349bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 68449bd79ccSDebojyoti Ghosh jrow = ii[i]; 68549bd79ccSDebojyoti Ghosh n = ii[i+1] - jrow; 68649bd79ccSDebojyoti Ghosh sums = y + p*i; 68749bd79ccSDebojyoti Ghosh for (j=0; j<n; j++) { 68849bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 68949bd79ccSDebojyoti Ghosh sums[k] += v[jrow+j]*x[q*idx[jrow+j]+k]; 69049bd79ccSDebojyoti Ghosh } 69149bd79ccSDebojyoti Ghosh } 69249bd79ccSDebojyoti Ghosh } 693ca0c957dSBarry Smith ierr = PetscLogFlops(3.0*(a->nz)*p);CHKERRQ(ierr); 69449bd79ccSDebojyoti Ghosh } else if (t) { 69549bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 69649bd79ccSDebojyoti Ghosh jrow = ii[i]; 69749bd79ccSDebojyoti Ghosh n = ii[i+1] - jrow; 69849bd79ccSDebojyoti Ghosh sums = y + p*i; 69949bd79ccSDebojyoti Ghosh for (j=0; j<n; j++) { 70049bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 70149bd79ccSDebojyoti Ghosh for (l=0; l<q; l++) { 70249bd79ccSDebojyoti Ghosh sums[k] += v[jrow+j]*t[k+l*p]*x[q*idx[jrow+j]+l]; 70349bd79ccSDebojyoti Ghosh } 70449bd79ccSDebojyoti Ghosh } 70549bd79ccSDebojyoti Ghosh } 70649bd79ccSDebojyoti Ghosh } 707234d9204SRichard Tran Mills /* The flop count below assumes that v[jrow+j] is hoisted out (which an optimizing compiler is likely to do), 708234d9204SRichard Tran Mills * and also that T part is hoisted outside this loop (in exchange for temporary storage) as (A \otimes I) (I \otimes T), 709234d9204SRichard Tran Mills * so that this multiply doesn't have to be redone for each matrix entry, but just once per column. The latter 710234d9204SRichard Tran Mills * transformation is much less likely to be applied, but we nonetheless count the minimum flops required. */ 711ca0c957dSBarry Smith ierr = PetscLogFlops((2.0*p*q-p)*m+2.0*p*a->nz);CHKERRQ(ierr); 71249bd79ccSDebojyoti Ghosh } 71349bd79ccSDebojyoti Ghosh if (s) { 71449bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 71549bd79ccSDebojyoti Ghosh sums = y + p*i; 71649bd79ccSDebojyoti Ghosh bx = x + q*i; 71749bd79ccSDebojyoti Ghosh if (i < b->AIJ->cmap->n) { 71849bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 71949bd79ccSDebojyoti Ghosh for (k=0; k<p; k++) { 72049bd79ccSDebojyoti Ghosh sums[k] += s[k+j*p]*bx[j]; 72149bd79ccSDebojyoti Ghosh } 72249bd79ccSDebojyoti Ghosh } 72349bd79ccSDebojyoti Ghosh } 72449bd79ccSDebojyoti Ghosh } 725ca0c957dSBarry Smith ierr = PetscLogFlops(2.0*m*p*q);CHKERRQ(ierr); 72649bd79ccSDebojyoti Ghosh } 72749bd79ccSDebojyoti Ghosh 72849bd79ccSDebojyoti Ghosh ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 72949bd79ccSDebojyoti Ghosh ierr = VecRestoreArray(zz,&y);CHKERRQ(ierr); 73049bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 73149bd79ccSDebojyoti Ghosh } 73249bd79ccSDebojyoti Ghosh 733bb6fb833SRichard Tran Mills PetscErrorCode MatMult_SeqKAIJ(Mat A,Vec xx,Vec yy) 73449bd79ccSDebojyoti Ghosh { 73549bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 73649bd79ccSDebojyoti Ghosh PetscFunctionBegin; 737836168d5SRichard Tran Mills ierr = MatMultAdd_SeqKAIJ(A,xx,PETSC_NULL,yy);CHKERRQ(ierr); 73849bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 73949bd79ccSDebojyoti Ghosh } 74049bd79ccSDebojyoti Ghosh 74149bd79ccSDebojyoti Ghosh #include <petsc/private/kernels/blockinvert.h> 74249bd79ccSDebojyoti Ghosh 743bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_SeqKAIJ(Mat A,const PetscScalar **values) 74449bd79ccSDebojyoti Ghosh { 74549bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*)A->data; 74649bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ*)b->AIJ->data; 74749bd79ccSDebojyoti Ghosh const PetscScalar *S = b->S; 74849bd79ccSDebojyoti Ghosh const PetscScalar *T = b->T; 74949bd79ccSDebojyoti Ghosh const PetscScalar *v = a->a; 75049bd79ccSDebojyoti Ghosh const PetscInt p = b->p, q = b->q, m = b->AIJ->rmap->n, *idx = a->j, *ii = a->i; 75149bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 75249bd79ccSDebojyoti Ghosh PetscInt i,j,*v_pivots,dof,dof2; 75349bd79ccSDebojyoti Ghosh PetscScalar *diag,aval,*v_work; 75449bd79ccSDebojyoti Ghosh 75549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 75649bd79ccSDebojyoti Ghosh if (p != q) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATKAIJ: Block size must be square to calculate inverse."); 75731a97b9aSRichard Tran Mills if ((!S) && (!T) && (!b->isTI)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATKAIJ: Cannot invert a zero matrix."); 75849bd79ccSDebojyoti Ghosh 75949bd79ccSDebojyoti Ghosh dof = p; 76049bd79ccSDebojyoti Ghosh dof2 = dof*dof; 76149bd79ccSDebojyoti Ghosh 76249bd79ccSDebojyoti Ghosh if (b->ibdiagvalid) { 76349bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 76449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 76549bd79ccSDebojyoti Ghosh } 76649bd79ccSDebojyoti Ghosh if (!b->ibdiag) { 767a84f8069SRichard Tran Mills ierr = PetscMalloc1(dof2*m*sizeof(PetscScalar),&b->ibdiag);CHKERRQ(ierr); 76849bd79ccSDebojyoti Ghosh ierr = PetscLogObjectMemory((PetscObject)A,dof2*m*sizeof(PetscScalar));CHKERRQ(ierr); 76949bd79ccSDebojyoti Ghosh } 77049bd79ccSDebojyoti Ghosh if (values) *values = b->ibdiag; 77149bd79ccSDebojyoti Ghosh diag = b->ibdiag; 77249bd79ccSDebojyoti Ghosh 77349bd79ccSDebojyoti Ghosh ierr = PetscMalloc2(dof,&v_work,dof,&v_pivots);CHKERRQ(ierr); 77449bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 77549bd79ccSDebojyoti Ghosh if (S) { 77649bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(diag,S,dof2*sizeof(PetscScalar));CHKERRQ(ierr); 77749bd79ccSDebojyoti Ghosh } else { 77849bd79ccSDebojyoti Ghosh ierr = PetscMemzero(diag,dof2*sizeof(PetscScalar));CHKERRQ(ierr); 77949bd79ccSDebojyoti Ghosh } 78049bd79ccSDebojyoti Ghosh if (b->isTI) { 78149bd79ccSDebojyoti Ghosh aval = 0; 78249bd79ccSDebojyoti Ghosh for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j]; 78349bd79ccSDebojyoti Ghosh for (j=0; j<dof; j++) diag[j+dof*j] += aval; 78449bd79ccSDebojyoti Ghosh } else if (T) { 78549bd79ccSDebojyoti Ghosh aval = 0; 78649bd79ccSDebojyoti Ghosh for (j=ii[i]; j<ii[i+1]; j++) if (idx[j] == i) aval = v[j]; 78749bd79ccSDebojyoti Ghosh for (j=0; j<dof2; j++) diag[j] += aval*T[j]; 78849bd79ccSDebojyoti Ghosh } 78949bd79ccSDebojyoti Ghosh ierr = PetscKernel_A_gets_inverse_A(dof,diag,v_pivots,v_work,PETSC_FALSE,NULL);CHKERRQ(ierr); 79049bd79ccSDebojyoti Ghosh diag += dof2; 79149bd79ccSDebojyoti Ghosh } 79249bd79ccSDebojyoti Ghosh ierr = PetscFree2(v_work,v_pivots);CHKERRQ(ierr); 79349bd79ccSDebojyoti Ghosh 79449bd79ccSDebojyoti Ghosh b->ibdiagvalid = PETSC_TRUE; 79549bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 79649bd79ccSDebojyoti Ghosh } 79749bd79ccSDebojyoti Ghosh 79849bd79ccSDebojyoti Ghosh static PetscErrorCode MatGetDiagonalBlock_MPIKAIJ(Mat A,Mat *B) 79949bd79ccSDebojyoti Ghosh { 80049bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *kaij = (Mat_MPIKAIJ*) A->data; 80149bd79ccSDebojyoti Ghosh 80249bd79ccSDebojyoti Ghosh PetscFunctionBegin; 80349bd79ccSDebojyoti Ghosh *B = kaij->AIJ; 80449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 80549bd79ccSDebojyoti Ghosh } 80649bd79ccSDebojyoti Ghosh 80749bd79ccSDebojyoti Ghosh PetscErrorCode MatSOR_SeqKAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 80849bd79ccSDebojyoti Ghosh { 80949bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 81049bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *kaij = (Mat_SeqKAIJ*) A->data; 81149bd79ccSDebojyoti Ghosh Mat_SeqAIJ *a = (Mat_SeqAIJ*)kaij->AIJ->data; 81249bd79ccSDebojyoti Ghosh const PetscScalar *aa = a->a, *T = kaij->T, *v; 81349bd79ccSDebojyoti Ghosh const PetscInt m = kaij->AIJ->rmap->n, *ai=a->i, *aj=a->j, p = kaij->p, q = kaij->q, *diag, *vi; 81449bd79ccSDebojyoti Ghosh const PetscScalar *b, *xb, *idiag; 81549bd79ccSDebojyoti Ghosh PetscScalar *x, *work, *workt, *w, *y, *arr, *t, *arrt; 81649bd79ccSDebojyoti Ghosh PetscInt i, j, k, i2, bs, bs2, nz; 81749bd79ccSDebojyoti Ghosh 81849bd79ccSDebojyoti Ghosh PetscFunctionBegin; 81949bd79ccSDebojyoti Ghosh its = its*lits; 82049bd79ccSDebojyoti Ghosh if (flag & SOR_EISENSTAT) SETERRQ (PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 821c0aa6a63SJacob Faibussowitsch if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 8226a375485SRichard Tran Mills if (fshift) SETERRQ (PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for diagonal shift"); 8236a375485SRichard Tran Mills if ((flag & SOR_APPLY_UPPER) || (flag & SOR_APPLY_LOWER)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for applying upper or lower triangular parts"); 8246a375485SRichard Tran Mills if (p != q) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSOR for KAIJ: No support for non-square dense blocks"); 82549bd79ccSDebojyoti Ghosh else {bs = p; bs2 = bs*bs; } 82649bd79ccSDebojyoti Ghosh 82749bd79ccSDebojyoti Ghosh if (!m) PetscFunctionReturn(0); 82849bd79ccSDebojyoti Ghosh 829bb6fb833SRichard Tran Mills if (!kaij->ibdiagvalid) { ierr = MatInvertBlockDiagonal_SeqKAIJ(A,NULL);CHKERRQ(ierr); } 83049bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 83149bd79ccSDebojyoti Ghosh diag = a->diag; 83249bd79ccSDebojyoti Ghosh 83349bd79ccSDebojyoti Ghosh if (!kaij->sor.setup) { 83449bd79ccSDebojyoti Ghosh ierr = PetscMalloc5(bs,&kaij->sor.w,bs,&kaij->sor.y,m*bs,&kaij->sor.work,m*bs,&kaij->sor.t,m*bs2,&kaij->sor.arr);CHKERRQ(ierr); 83549bd79ccSDebojyoti Ghosh kaij->sor.setup = PETSC_TRUE; 83649bd79ccSDebojyoti Ghosh } 83749bd79ccSDebojyoti Ghosh y = kaij->sor.y; 83849bd79ccSDebojyoti Ghosh w = kaij->sor.w; 83949bd79ccSDebojyoti Ghosh work = kaij->sor.work; 84049bd79ccSDebojyoti Ghosh t = kaij->sor.t; 84149bd79ccSDebojyoti Ghosh arr = kaij->sor.arr; 84249bd79ccSDebojyoti Ghosh 84349bd79ccSDebojyoti Ghosh ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 84449bd79ccSDebojyoti Ghosh ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 84549bd79ccSDebojyoti Ghosh 84649bd79ccSDebojyoti Ghosh if (flag & SOR_ZERO_INITIAL_GUESS) { 84749bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 84849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,b,idiag,x); /* x[0:bs] <- D^{-1} b[0:bs] */ 84949bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); 85049bd79ccSDebojyoti Ghosh i2 = bs; 85149bd79ccSDebojyoti Ghosh idiag += bs2; 85249bd79ccSDebojyoti Ghosh for (i=1; i<m; i++) { 85349bd79ccSDebojyoti Ghosh v = aa + ai[i]; 85449bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 85549bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 85649bd79ccSDebojyoti Ghosh 85749bd79ccSDebojyoti Ghosh if (T) { /* b - T (Arow * x) */ 8582ae760e3SRichard Tran Mills ierr = PetscMemzero(w,bs*sizeof(PetscScalar));CHKERRQ(ierr); 85949bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 86049bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k]; 86149bd79ccSDebojyoti Ghosh } 86249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs,w,T,&t[i2]); 86349bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) t[i2+k] += b[i2+k]; 86449bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 8659eb573c1SRichard Tran Mills ierr = PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 86649bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 86749bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) t[i2+k] -= v[j] * x[vi[j]*bs+k]; 86849bd79ccSDebojyoti Ghosh } 86949bd79ccSDebojyoti Ghosh } else { 8709eb573c1SRichard Tran Mills ierr = PetscMemcpy(t+i2,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 87149bd79ccSDebojyoti Ghosh } 87249bd79ccSDebojyoti Ghosh 87349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,t+i2,idiag,y); 87449bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) x[i2+j] = omega * y[j]; 87549bd79ccSDebojyoti Ghosh 87649bd79ccSDebojyoti Ghosh idiag += bs2; 87749bd79ccSDebojyoti Ghosh i2 += bs; 87849bd79ccSDebojyoti Ghosh } 87949bd79ccSDebojyoti Ghosh /* for logging purposes assume number of nonzero in lower half is 1/2 of total */ 88049bd79ccSDebojyoti Ghosh ierr = PetscLogFlops(1.0*bs2*a->nz);CHKERRQ(ierr); 88149bd79ccSDebojyoti Ghosh xb = t; 88249bd79ccSDebojyoti Ghosh } else xb = b; 88349bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 88449bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag+bs2*(m-1); 88549bd79ccSDebojyoti Ghosh i2 = bs * (m-1); 88649bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 88749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,x+i2); 88849bd79ccSDebojyoti Ghosh i2 -= bs; 88949bd79ccSDebojyoti Ghosh idiag -= bs2; 89049bd79ccSDebojyoti Ghosh for (i=m-2; i>=0; i--) { 89149bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1 ; 89249bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 89349bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 89449bd79ccSDebojyoti Ghosh 89549bd79ccSDebojyoti Ghosh if (T) { /* FIXME: This branch untested */ 89649bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 89749bd79ccSDebojyoti Ghosh /* copy all rows of x that are needed into contiguous space */ 89849bd79ccSDebojyoti Ghosh workt = work; 89949bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 90049bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 90149bd79ccSDebojyoti Ghosh workt += bs; 90249bd79ccSDebojyoti Ghosh } 90349bd79ccSDebojyoti Ghosh arrt = arr; 90449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 90549bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 90649bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 90749bd79ccSDebojyoti Ghosh arrt += bs2; 90849bd79ccSDebojyoti Ghosh } 90949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 91049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 9119eb573c1SRichard Tran Mills ierr = PetscMemcpy(w,t+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 91249bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 91349bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) w[k] -= v[j] * x[vi[j]*bs+k]; 91449bd79ccSDebojyoti Ghosh } 91549bd79ccSDebojyoti Ghosh } 91649bd79ccSDebojyoti Ghosh 91749bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); /* RHS incorrect for omega != 1.0 */ 91849bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) x[i2+j] = (1.0-omega) * x[i2+j] + omega * y[j]; 91949bd79ccSDebojyoti Ghosh 92049bd79ccSDebojyoti Ghosh idiag -= bs2; 92149bd79ccSDebojyoti Ghosh i2 -= bs; 92249bd79ccSDebojyoti Ghosh } 92349bd79ccSDebojyoti Ghosh ierr = PetscLogFlops(1.0*bs2*(a->nz));CHKERRQ(ierr); 92449bd79ccSDebojyoti Ghosh } 92549bd79ccSDebojyoti Ghosh its--; 92649bd79ccSDebojyoti Ghosh } 92749bd79ccSDebojyoti Ghosh while (its--) { /* FIXME: This branch not updated */ 92849bd79ccSDebojyoti Ghosh if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 92949bd79ccSDebojyoti Ghosh i2 = 0; 93049bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag; 93149bd79ccSDebojyoti Ghosh for (i=0; i<m; i++) { 93249bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 93349bd79ccSDebojyoti Ghosh 93449bd79ccSDebojyoti Ghosh v = aa + ai[i]; 93549bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 93649bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 93749bd79ccSDebojyoti Ghosh workt = work; 93849bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 93949bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 94049bd79ccSDebojyoti Ghosh workt += bs; 94149bd79ccSDebojyoti Ghosh } 94249bd79ccSDebojyoti Ghosh arrt = arr; 94349bd79ccSDebojyoti Ghosh if (T) { 94449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 94549bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 94649bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 94749bd79ccSDebojyoti Ghosh arrt += bs2; 94849bd79ccSDebojyoti Ghosh } 94949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 95049bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 95149bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 95249bd79ccSDebojyoti Ghosh ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 95349bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 95449bd79ccSDebojyoti Ghosh arrt += bs2; 95549bd79ccSDebojyoti Ghosh } 95649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 95749bd79ccSDebojyoti Ghosh } 95849bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(t+i2,w,bs*sizeof(PetscScalar));CHKERRQ(ierr); 95949bd79ccSDebojyoti Ghosh 96049bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 96149bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 96249bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 96349bd79ccSDebojyoti Ghosh workt = work; 96449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 96549bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 96649bd79ccSDebojyoti Ghosh workt += bs; 96749bd79ccSDebojyoti Ghosh } 96849bd79ccSDebojyoti Ghosh arrt = arr; 96949bd79ccSDebojyoti Ghosh if (T) { 97049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 97149bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 97249bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 97349bd79ccSDebojyoti Ghosh arrt += bs2; 97449bd79ccSDebojyoti Ghosh } 97549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 97649bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 97749bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 97849bd79ccSDebojyoti Ghosh ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 97949bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 98049bd79ccSDebojyoti Ghosh arrt += bs2; 98149bd79ccSDebojyoti Ghosh } 98249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 98349bd79ccSDebojyoti Ghosh } 98449bd79ccSDebojyoti Ghosh 98549bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 98649bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 98749bd79ccSDebojyoti Ghosh 98849bd79ccSDebojyoti Ghosh idiag += bs2; 98949bd79ccSDebojyoti Ghosh i2 += bs; 99049bd79ccSDebojyoti Ghosh } 99149bd79ccSDebojyoti Ghosh xb = t; 99249bd79ccSDebojyoti Ghosh } 99349bd79ccSDebojyoti Ghosh else xb = b; 99449bd79ccSDebojyoti Ghosh if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 99549bd79ccSDebojyoti Ghosh idiag = kaij->ibdiag+bs2*(m-1); 99649bd79ccSDebojyoti Ghosh i2 = bs * (m-1); 99749bd79ccSDebojyoti Ghosh if (xb == b) { 99849bd79ccSDebojyoti Ghosh for (i=m-1; i>=0; i--) { 99949bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(w,b+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 100049bd79ccSDebojyoti Ghosh 100149bd79ccSDebojyoti Ghosh v = aa + ai[i]; 100249bd79ccSDebojyoti Ghosh vi = aj + ai[i]; 100349bd79ccSDebojyoti Ghosh nz = diag[i] - ai[i]; 100449bd79ccSDebojyoti Ghosh workt = work; 100549bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 100649bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 100749bd79ccSDebojyoti Ghosh workt += bs; 100849bd79ccSDebojyoti Ghosh } 100949bd79ccSDebojyoti Ghosh arrt = arr; 101049bd79ccSDebojyoti Ghosh if (T) { 101149bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 101249bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 101349bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 101449bd79ccSDebojyoti Ghosh arrt += bs2; 101549bd79ccSDebojyoti Ghosh } 101649bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 101749bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 101849bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 101949bd79ccSDebojyoti Ghosh ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 102049bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 102149bd79ccSDebojyoti Ghosh arrt += bs2; 102249bd79ccSDebojyoti Ghosh } 102349bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 102449bd79ccSDebojyoti Ghosh } 102549bd79ccSDebojyoti Ghosh 102649bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 102749bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 102849bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 102949bd79ccSDebojyoti Ghosh workt = work; 103049bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 103149bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 103249bd79ccSDebojyoti Ghosh workt += bs; 103349bd79ccSDebojyoti Ghosh } 103449bd79ccSDebojyoti Ghosh arrt = arr; 103549bd79ccSDebojyoti Ghosh if (T) { 103649bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 103749bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 103849bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 103949bd79ccSDebojyoti Ghosh arrt += bs2; 104049bd79ccSDebojyoti Ghosh } 104149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 104249bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 104349bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 104449bd79ccSDebojyoti Ghosh ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 104549bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 104649bd79ccSDebojyoti Ghosh arrt += bs2; 104749bd79ccSDebojyoti Ghosh } 104849bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 104949bd79ccSDebojyoti Ghosh } 105049bd79ccSDebojyoti Ghosh 105149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 105249bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 105349bd79ccSDebojyoti Ghosh } 105449bd79ccSDebojyoti Ghosh } else { 105549bd79ccSDebojyoti Ghosh for (i=m-1; i>=0; i--) { 105649bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(w,xb+i2,bs*sizeof(PetscScalar));CHKERRQ(ierr); 105749bd79ccSDebojyoti Ghosh v = aa + diag[i] + 1; 105849bd79ccSDebojyoti Ghosh vi = aj + diag[i] + 1; 105949bd79ccSDebojyoti Ghosh nz = ai[i+1] - diag[i] - 1; 106049bd79ccSDebojyoti Ghosh workt = work; 106149bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 106249bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(workt,x + bs*(*vi++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 106349bd79ccSDebojyoti Ghosh workt += bs; 106449bd79ccSDebojyoti Ghosh } 106549bd79ccSDebojyoti Ghosh arrt = arr; 106649bd79ccSDebojyoti Ghosh if (T) { 106749bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 106849bd79ccSDebojyoti Ghosh ierr = PetscMemcpy(arrt,T,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 106949bd79ccSDebojyoti Ghosh for (k=0; k<bs2; k++) arrt[k] *= v[j]; 107049bd79ccSDebojyoti Ghosh arrt += bs2; 107149bd79ccSDebojyoti Ghosh } 107249bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 107349bd79ccSDebojyoti Ghosh } else if (kaij->isTI) { 107449bd79ccSDebojyoti Ghosh for (j=0; j<nz; j++) { 107549bd79ccSDebojyoti Ghosh ierr = PetscMemzero(arrt,bs2*sizeof(PetscScalar));CHKERRQ(ierr); 107649bd79ccSDebojyoti Ghosh for (k=0; k<bs; k++) arrt[k+bs*k] = v[j]; 107749bd79ccSDebojyoti Ghosh arrt += bs2; 107849bd79ccSDebojyoti Ghosh } 107949bd79ccSDebojyoti Ghosh PetscKernel_w_gets_w_minus_Ar_times_v(bs,bs*nz,w,arr,work); 108049bd79ccSDebojyoti Ghosh } 108149bd79ccSDebojyoti Ghosh PetscKernel_w_gets_Ar_times_v(bs,bs,w,idiag,y); 108249bd79ccSDebojyoti Ghosh for (j=0; j<bs; j++) *(x+i2+j) = (1.0-omega) * *(x+i2+j) + omega * *(y+j); 108349bd79ccSDebojyoti Ghosh } 108449bd79ccSDebojyoti Ghosh } 108549bd79ccSDebojyoti Ghosh ierr = PetscLogFlops(1.0*bs2*(a->nz));CHKERRQ(ierr); 108649bd79ccSDebojyoti Ghosh } 108749bd79ccSDebojyoti Ghosh } 108849bd79ccSDebojyoti Ghosh 108949bd79ccSDebojyoti Ghosh ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 109049bd79ccSDebojyoti Ghosh ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 109149bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 109249bd79ccSDebojyoti Ghosh } 109349bd79ccSDebojyoti Ghosh 109449bd79ccSDebojyoti Ghosh /*===================================================================================*/ 109549bd79ccSDebojyoti Ghosh 1096836168d5SRichard Tran Mills PetscErrorCode MatMultAdd_MPIKAIJ(Mat A,Vec xx,Vec yy,Vec zz) 109749bd79ccSDebojyoti Ghosh { 109849bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 109949bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 110049bd79ccSDebojyoti Ghosh 110149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 110249bd79ccSDebojyoti Ghosh if (!yy) { 110349bd79ccSDebojyoti Ghosh ierr = VecSet(zz,0.0);CHKERRQ(ierr); 110449bd79ccSDebojyoti Ghosh } else { 110549bd79ccSDebojyoti Ghosh ierr = VecCopy(yy,zz);CHKERRQ(ierr); 110649bd79ccSDebojyoti Ghosh } 1107761d359dSRichard Tran Mills ierr = MatKAIJ_build_AIJ_OAIJ(A);CHKERRQ(ierr); /* Ensure b->AIJ and b->OAIJ are up to date. */ 110849bd79ccSDebojyoti Ghosh /* start the scatter */ 110949bd79ccSDebojyoti Ghosh ierr = VecScatterBegin(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 111049bd79ccSDebojyoti Ghosh ierr = (*b->AIJ->ops->multadd)(b->AIJ,xx,zz,zz);CHKERRQ(ierr); 111149bd79ccSDebojyoti Ghosh ierr = VecScatterEnd(b->ctx,xx,b->w,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 111249bd79ccSDebojyoti Ghosh ierr = (*b->OAIJ->ops->multadd)(b->OAIJ,b->w,zz,zz);CHKERRQ(ierr); 111349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 111449bd79ccSDebojyoti Ghosh } 111549bd79ccSDebojyoti Ghosh 1116bb6fb833SRichard Tran Mills PetscErrorCode MatMult_MPIKAIJ(Mat A,Vec xx,Vec yy) 111749bd79ccSDebojyoti Ghosh { 111849bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 111949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1120836168d5SRichard Tran Mills ierr = MatMultAdd_MPIKAIJ(A,xx,PETSC_NULL,yy);CHKERRQ(ierr); 112149bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 112249bd79ccSDebojyoti Ghosh } 112349bd79ccSDebojyoti Ghosh 1124bb6fb833SRichard Tran Mills PetscErrorCode MatInvertBlockDiagonal_MPIKAIJ(Mat A,const PetscScalar **values) 112549bd79ccSDebojyoti Ghosh { 112649bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*)A->data; 112749bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 112849bd79ccSDebojyoti Ghosh 112949bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1130761d359dSRichard Tran Mills ierr = MatKAIJ_build_AIJ_OAIJ(A);CHKERRQ(ierr); /* Ensure b->AIJ is up to date. */ 113149bd79ccSDebojyoti Ghosh ierr = (*b->AIJ->ops->invertblockdiagonal)(b->AIJ,values);CHKERRQ(ierr); 113249bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 113349bd79ccSDebojyoti Ghosh } 113449bd79ccSDebojyoti Ghosh 113549bd79ccSDebojyoti Ghosh /* ----------------------------------------------------------------*/ 113649bd79ccSDebojyoti Ghosh 113749bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values) 113849bd79ccSDebojyoti Ghosh { 113949bd79ccSDebojyoti Ghosh Mat_SeqKAIJ *b = (Mat_SeqKAIJ*) A->data; 11401ca5ffdbSRichard Tran Mills PetscErrorCode diag = PETSC_FALSE; 11411ca5ffdbSRichard Tran Mills PetscErrorCode ierr; 114249bd79ccSDebojyoti Ghosh PetscInt nzaij,nz,*colsaij,*idx,i,j,p=b->p,q=b->q,r=row/p,s=row%p,c; 114349bd79ccSDebojyoti Ghosh PetscScalar *vaij,*v,*S=b->S,*T=b->T; 114449bd79ccSDebojyoti Ghosh 114549bd79ccSDebojyoti Ghosh PetscFunctionBegin; 114649bd79ccSDebojyoti Ghosh if (b->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active"); 114749bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 1148c0aa6a63SJacob Faibussowitsch if (row < 0 || row >= A->rmap->n) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " out of range",row); 114949bd79ccSDebojyoti Ghosh 115049bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 115149bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 115249bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 115349bd79ccSDebojyoti Ghosh if (values) *values = NULL; 115449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 115549bd79ccSDebojyoti Ghosh } 115649bd79ccSDebojyoti Ghosh 115749bd79ccSDebojyoti Ghosh if (T || b->isTI) { 115849bd79ccSDebojyoti Ghosh ierr = MatGetRow_SeqAIJ(b->AIJ,r,&nzaij,&colsaij,&vaij);CHKERRQ(ierr); 115949bd79ccSDebojyoti Ghosh c = nzaij; 116049bd79ccSDebojyoti Ghosh for (i=0; i<nzaij; i++) { 116149bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 116249bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 116349bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 116449bd79ccSDebojyoti Ghosh c = i; 116549bd79ccSDebojyoti Ghosh } 116649bd79ccSDebojyoti Ghosh } 116749bd79ccSDebojyoti Ghosh } else nzaij = c = 0; 116849bd79ccSDebojyoti Ghosh 116949bd79ccSDebojyoti Ghosh /* calculate size of row */ 117049bd79ccSDebojyoti Ghosh nz = 0; 117149bd79ccSDebojyoti Ghosh if (S) nz += q; 117249bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (nzaij-1)*q : nzaij*q); 117349bd79ccSDebojyoti Ghosh 117449bd79ccSDebojyoti Ghosh if (cols || values) { 117549bd79ccSDebojyoti Ghosh ierr = PetscMalloc2(nz,&idx,nz,&v);CHKERRQ(ierr); 117638822f9dSRichard Tran Mills for (i=0; i<q; i++) { 117738822f9dSRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 117838822f9dSRichard Tran Mills v[i] = 0.0; 117938822f9dSRichard Tran Mills } 118049bd79ccSDebojyoti Ghosh if (b->isTI) { 118149bd79ccSDebojyoti Ghosh for (i=0; i<nzaij; i++) { 118249bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 118349bd79ccSDebojyoti Ghosh idx[i*q+j] = colsaij[i]*q+j; 118449bd79ccSDebojyoti Ghosh v[i*q+j] = (j==s ? vaij[i] : 0); 118549bd79ccSDebojyoti Ghosh } 118649bd79ccSDebojyoti Ghosh } 118749bd79ccSDebojyoti Ghosh } else if (T) { 118849bd79ccSDebojyoti Ghosh for (i=0; i<nzaij; i++) { 118949bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 119049bd79ccSDebojyoti Ghosh idx[i*q+j] = colsaij[i]*q+j; 119149bd79ccSDebojyoti Ghosh v[i*q+j] = vaij[i]*T[s+j*p]; 119249bd79ccSDebojyoti Ghosh } 119349bd79ccSDebojyoti Ghosh } 119449bd79ccSDebojyoti Ghosh } 119549bd79ccSDebojyoti Ghosh if (S) { 119649bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 119749bd79ccSDebojyoti Ghosh idx[c*q+j] = r*q+j; 119849bd79ccSDebojyoti Ghosh v[c*q+j] += S[s+j*p]; 119949bd79ccSDebojyoti Ghosh } 120049bd79ccSDebojyoti Ghosh } 120149bd79ccSDebojyoti Ghosh } 120249bd79ccSDebojyoti Ghosh 120349bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 120449bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 120549bd79ccSDebojyoti Ghosh if (values) *values = v; 120649bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 120749bd79ccSDebojyoti Ghosh } 120849bd79ccSDebojyoti Ghosh 120949bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_SeqKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 121049bd79ccSDebojyoti Ghosh { 121149bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 1212cb4a9cd9SHong Zhang 121349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1214cb4a9cd9SHong Zhang if (nz) *nz = 0; 121549bd79ccSDebojyoti Ghosh ierr = PetscFree2(*idx,*v);CHKERRQ(ierr); 121649bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE; 121749bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 121849bd79ccSDebojyoti Ghosh } 121949bd79ccSDebojyoti Ghosh 122049bd79ccSDebojyoti Ghosh PetscErrorCode MatGetRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **values) 122149bd79ccSDebojyoti Ghosh { 122249bd79ccSDebojyoti Ghosh Mat_MPIKAIJ *b = (Mat_MPIKAIJ*) A->data; 122349bd79ccSDebojyoti Ghosh Mat AIJ = b->A; 1224fc64b2cfSRichard Tran Mills PetscBool diag = PETSC_FALSE; 1225761d359dSRichard Tran Mills Mat MatAIJ,MatOAIJ; 122649bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 122749bd79ccSDebojyoti Ghosh const PetscInt rstart=A->rmap->rstart,rend=A->rmap->rend,p=b->p,q=b->q,*garray; 1228389eba51SJed Brown PetscInt nz,*idx,ncolsaij = 0,ncolsoaij = 0,*colsaij,*colsoaij,r,s,c,i,j,lrow; 122949bd79ccSDebojyoti Ghosh PetscScalar *v,*vals,*ovals,*S=b->S,*T=b->T; 123049bd79ccSDebojyoti Ghosh 123149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 1232761d359dSRichard Tran Mills ierr = MatKAIJ_build_AIJ_OAIJ(A);CHKERRQ(ierr); /* Ensure b->AIJ and b->OAIJ are up to date. */ 1233761d359dSRichard Tran Mills MatAIJ = ((Mat_SeqKAIJ*)b->AIJ->data)->AIJ; 1234761d359dSRichard Tran Mills MatOAIJ = ((Mat_SeqKAIJ*)b->OAIJ->data)->AIJ; 123549bd79ccSDebojyoti Ghosh if (b->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active"); 123649bd79ccSDebojyoti Ghosh b->getrowactive = PETSC_TRUE; 123749bd79ccSDebojyoti Ghosh if (row < rstart || row >= rend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only local rows"); 123849bd79ccSDebojyoti Ghosh lrow = row - rstart; 123949bd79ccSDebojyoti Ghosh 124049bd79ccSDebojyoti Ghosh if ((!S) && (!T) && (!b->isTI)) { 124149bd79ccSDebojyoti Ghosh if (ncols) *ncols = 0; 124249bd79ccSDebojyoti Ghosh if (cols) *cols = NULL; 124349bd79ccSDebojyoti Ghosh if (values) *values = NULL; 124449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 124549bd79ccSDebojyoti Ghosh } 124649bd79ccSDebojyoti Ghosh 124749bd79ccSDebojyoti Ghosh r = lrow/p; 124849bd79ccSDebojyoti Ghosh s = lrow%p; 124949bd79ccSDebojyoti Ghosh 125049bd79ccSDebojyoti Ghosh if (T || b->isTI) { 12512ae760e3SRichard Tran Mills ierr = MatMPIAIJGetSeqAIJ(AIJ,NULL,NULL,&garray);CHKERRQ(ierr); 125249bd79ccSDebojyoti Ghosh ierr = MatGetRow_SeqAIJ(MatAIJ,lrow/p,&ncolsaij,&colsaij,&vals);CHKERRQ(ierr); 125349bd79ccSDebojyoti Ghosh ierr = MatGetRow_SeqAIJ(MatOAIJ,lrow/p,&ncolsoaij,&colsoaij,&ovals);CHKERRQ(ierr); 125449bd79ccSDebojyoti Ghosh 125549bd79ccSDebojyoti Ghosh c = ncolsaij + ncolsoaij; 125649bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 125749bd79ccSDebojyoti Ghosh /* check if this row contains a diagonal entry */ 125849bd79ccSDebojyoti Ghosh if (colsaij[i] == r) { 125949bd79ccSDebojyoti Ghosh diag = PETSC_TRUE; 126049bd79ccSDebojyoti Ghosh c = i; 126149bd79ccSDebojyoti Ghosh } 126249bd79ccSDebojyoti Ghosh } 126349bd79ccSDebojyoti Ghosh } else c = 0; 126449bd79ccSDebojyoti Ghosh 126549bd79ccSDebojyoti Ghosh /* calculate size of row */ 126649bd79ccSDebojyoti Ghosh nz = 0; 126749bd79ccSDebojyoti Ghosh if (S) nz += q; 126849bd79ccSDebojyoti Ghosh if (T || b->isTI) nz += (diag && S ? (ncolsaij+ncolsoaij-1)*q : (ncolsaij+ncolsoaij)*q); 126949bd79ccSDebojyoti Ghosh 127049bd79ccSDebojyoti Ghosh if (cols || values) { 127149bd79ccSDebojyoti Ghosh ierr = PetscMalloc2(nz,&idx,nz,&v);CHKERRQ(ierr); 1272a437a796SRichard Tran Mills for (i=0; i<q; i++) { 1273a437a796SRichard Tran Mills /* We need to initialize the v[i] to zero to handle the case in which T is NULL (not the identity matrix). */ 1274a437a796SRichard Tran Mills v[i] = 0.0; 1275a437a796SRichard Tran Mills } 127649bd79ccSDebojyoti Ghosh if (b->isTI) { 127749bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 127849bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 127949bd79ccSDebojyoti Ghosh idx[i*q+j] = (colsaij[i]+rstart/p)*q+j; 128049bd79ccSDebojyoti Ghosh v[i*q+j] = (j==s ? vals[i] : 0.0); 128149bd79ccSDebojyoti Ghosh } 128249bd79ccSDebojyoti Ghosh } 128349bd79ccSDebojyoti Ghosh for (i=0; i<ncolsoaij; i++) { 128449bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 128549bd79ccSDebojyoti Ghosh idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j; 128649bd79ccSDebojyoti Ghosh v[(i+ncolsaij)*q+j] = (j==s ? ovals[i]: 0.0); 128749bd79ccSDebojyoti Ghosh } 128849bd79ccSDebojyoti Ghosh } 128949bd79ccSDebojyoti Ghosh } else if (T) { 129049bd79ccSDebojyoti Ghosh for (i=0; i<ncolsaij; i++) { 129149bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 129249bd79ccSDebojyoti Ghosh idx[i*q+j] = (colsaij[i]+rstart/p)*q+j; 129349bd79ccSDebojyoti Ghosh v[i*q+j] = vals[i]*T[s+j*p]; 129449bd79ccSDebojyoti Ghosh } 129549bd79ccSDebojyoti Ghosh } 129649bd79ccSDebojyoti Ghosh for (i=0; i<ncolsoaij; i++) { 129749bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 129849bd79ccSDebojyoti Ghosh idx[(i+ncolsaij)*q+j] = garray[colsoaij[i]]*q+j; 129949bd79ccSDebojyoti Ghosh v[(i+ncolsaij)*q+j] = ovals[i]*T[s+j*p]; 130049bd79ccSDebojyoti Ghosh } 130149bd79ccSDebojyoti Ghosh } 130249bd79ccSDebojyoti Ghosh } 130349bd79ccSDebojyoti Ghosh if (S) { 130449bd79ccSDebojyoti Ghosh for (j=0; j<q; j++) { 130549bd79ccSDebojyoti Ghosh idx[c*q+j] = (r+rstart/p)*q+j; 130649bd79ccSDebojyoti Ghosh v[c*q+j] += S[s+j*p]; 130749bd79ccSDebojyoti Ghosh } 130849bd79ccSDebojyoti Ghosh } 130949bd79ccSDebojyoti Ghosh } 131049bd79ccSDebojyoti Ghosh 131149bd79ccSDebojyoti Ghosh if (ncols) *ncols = nz; 131249bd79ccSDebojyoti Ghosh if (cols) *cols = idx; 131349bd79ccSDebojyoti Ghosh if (values) *values = v; 131449bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 131549bd79ccSDebojyoti Ghosh } 131649bd79ccSDebojyoti Ghosh 131749bd79ccSDebojyoti Ghosh PetscErrorCode MatRestoreRow_MPIKAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 131849bd79ccSDebojyoti Ghosh { 131949bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 132049bd79ccSDebojyoti Ghosh PetscFunctionBegin; 132149bd79ccSDebojyoti Ghosh ierr = PetscFree2(*idx,*v);CHKERRQ(ierr); 132249bd79ccSDebojyoti Ghosh ((Mat_SeqKAIJ*)A->data)->getrowactive = PETSC_FALSE; 132349bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 132449bd79ccSDebojyoti Ghosh } 132549bd79ccSDebojyoti Ghosh 132649bd79ccSDebojyoti Ghosh PetscErrorCode MatCreateSubMatrix_KAIJ(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) 132749bd79ccSDebojyoti Ghosh { 132849bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 132949bd79ccSDebojyoti Ghosh Mat A; 133049bd79ccSDebojyoti Ghosh 133149bd79ccSDebojyoti Ghosh PetscFunctionBegin; 133249bd79ccSDebojyoti Ghosh ierr = MatConvert(mat,MATAIJ,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 133349bd79ccSDebojyoti Ghosh ierr = MatCreateSubMatrix(A,isrow,iscol,cll,newmat);CHKERRQ(ierr); 133449bd79ccSDebojyoti Ghosh ierr = MatDestroy(&A);CHKERRQ(ierr); 133549bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 133649bd79ccSDebojyoti Ghosh } 133749bd79ccSDebojyoti Ghosh 133849bd79ccSDebojyoti Ghosh /* ---------------------------------------------------------------------------------- */ 133949bd79ccSDebojyoti Ghosh /*@C 134049bd79ccSDebojyoti Ghosh MatCreateKAIJ - Creates a matrix type to be used for matrices of the following form: 134149bd79ccSDebojyoti Ghosh 134249bd79ccSDebojyoti Ghosh [I \otimes S + A \otimes T] 134349bd79ccSDebojyoti Ghosh 134449bd79ccSDebojyoti Ghosh where 134549bd79ccSDebojyoti Ghosh S is a dense (p \times q) matrix 134649bd79ccSDebojyoti Ghosh T is a dense (p \times q) matrix 134749bd79ccSDebojyoti Ghosh A is an AIJ (n \times n) matrix 134849bd79ccSDebojyoti Ghosh I is the identity matrix 134949bd79ccSDebojyoti Ghosh The resulting matrix is (np \times nq) 135049bd79ccSDebojyoti Ghosh 1351d3b90ce1SRichard Tran Mills S and T are always stored independently on all processes as PetscScalar arrays in column-major format. 135249bd79ccSDebojyoti Ghosh 135349bd79ccSDebojyoti Ghosh Collective 135449bd79ccSDebojyoti Ghosh 135549bd79ccSDebojyoti Ghosh Input Parameters: 135649bd79ccSDebojyoti Ghosh + A - the AIJ matrix 135749bd79ccSDebojyoti Ghosh . p - number of rows in S and T 1358d3b90ce1SRichard Tran Mills . q - number of columns in S and T 1359d3b90ce1SRichard Tran Mills . S - the S matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major) 1360d3b90ce1SRichard Tran Mills - T - the T matrix (can be PETSC_NULL), stored as a PetscScalar array (column-major) 136149bd79ccSDebojyoti Ghosh 136249bd79ccSDebojyoti Ghosh Output Parameter: 136349bd79ccSDebojyoti Ghosh . kaij - the new KAIJ matrix 136449bd79ccSDebojyoti Ghosh 1365d3b90ce1SRichard Tran Mills Notes: 1366d3b90ce1SRichard 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. 1367d3b90ce1SRichard Tran Mills Changes to the entries of the AIJ matrix will immediately affect the KAIJ matrix. 136849bd79ccSDebojyoti Ghosh 1369761d359dSRichard Tran Mills Developer Notes: 1370761d359dSRichard Tran Mills In the MATMPIKAIJ case, the internal 'AIJ' and 'OAIJ' sequential KAIJ matrices are kept up to date by tracking the object state 1371761d359dSRichard Tran Mills of the AIJ matrix 'A' that describes the blockwise action of the MATMPIKAIJ matrix and, if the object state has changed, lazily 1372761d359dSRichard Tran Mills rebuilding 'AIJ' and 'OAIJ' just before executing operations with the MATMPIKAIJ matrix. If new types of operations are added, 1373761d359dSRichard Tran Mills routines implementing those must also ensure these are rebuilt when needed (by calling the internal MatKAIJ_build_AIJ_OAIJ() routine). 1374761d359dSRichard Tran Mills 137549bd79ccSDebojyoti Ghosh Level: advanced 137649bd79ccSDebojyoti Ghosh 13770567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MATKAIJ 137849bd79ccSDebojyoti Ghosh @*/ 137949bd79ccSDebojyoti Ghosh PetscErrorCode MatCreateKAIJ(Mat A,PetscInt p,PetscInt q,const PetscScalar S[],const PetscScalar T[],Mat *kaij) 138049bd79ccSDebojyoti Ghosh { 138149bd79ccSDebojyoti Ghosh PetscErrorCode ierr; 138249bd79ccSDebojyoti Ghosh 138349bd79ccSDebojyoti Ghosh PetscFunctionBegin; 13840567c835SRichard Tran Mills ierr = MatCreate(PetscObjectComm((PetscObject)A),kaij);CHKERRQ(ierr); 1385*06d61a37SPierre Jolivet ierr = MatSetType(*kaij,MATKAIJ);CHKERRQ(ierr); 13860567c835SRichard Tran Mills ierr = MatKAIJSetAIJ(*kaij,A);CHKERRQ(ierr); 13870567c835SRichard Tran Mills ierr = MatKAIJSetS(*kaij,p,q,S);CHKERRQ(ierr); 13880567c835SRichard Tran Mills ierr = MatKAIJSetT(*kaij,p,q,T);CHKERRQ(ierr); 13892ae760e3SRichard Tran Mills ierr = MatSetUp(*kaij);CHKERRQ(ierr); 13900567c835SRichard Tran Mills PetscFunctionReturn(0); 13910567c835SRichard Tran Mills } 139249bd79ccSDebojyoti Ghosh 13930567c835SRichard Tran Mills /*MC 13945881e567SRichard Tran Mills MATKAIJ - MATKAIJ = "kaij" - A matrix type to be used to evaluate matrices of form 13955881e567SRichard Tran Mills [I \otimes S + A \otimes T], 13960567c835SRichard Tran Mills where 13975881e567SRichard Tran Mills S is a dense (p \times q) matrix, 13985881e567SRichard Tran Mills T is a dense (p \times q) matrix, 13995881e567SRichard Tran Mills A is an AIJ (n \times n) matrix, 14005881e567SRichard Tran Mills and I is the identity matrix. 14015881e567SRichard Tran Mills The resulting matrix is (np \times nq). 14020567c835SRichard Tran Mills 1403d3b90ce1SRichard Tran Mills S and T are always stored independently on all processes as PetscScalar arrays in column-major format. 14040567c835SRichard Tran Mills 14055881e567SRichard Tran Mills Notes: 14065881e567SRichard 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, 14075881e567SRichard Tran Mills where x and b are column vectors containing the row-major representations of X and B. 14085881e567SRichard Tran Mills 14090567c835SRichard Tran Mills Level: advanced 14100567c835SRichard Tran Mills 14110567c835SRichard Tran Mills .seealso: MatKAIJSetAIJ(), MatKAIJSetS(), MatKAIJSetT(), MatKAIJGetAIJ(), MatKAIJGetS(), MatKAIJGetT(), MatCreateKAIJ() 14120567c835SRichard Tran Mills M*/ 14130567c835SRichard Tran Mills 14140567c835SRichard Tran Mills PETSC_EXTERN PetscErrorCode MatCreate_KAIJ(Mat A) 14150567c835SRichard Tran Mills { 14160567c835SRichard Tran Mills PetscErrorCode ierr; 14170567c835SRichard Tran Mills Mat_MPIKAIJ *b; 14180567c835SRichard Tran Mills PetscMPIInt size; 14190567c835SRichard Tran Mills 14200567c835SRichard Tran Mills PetscFunctionBegin; 14210567c835SRichard Tran Mills ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 14220567c835SRichard Tran Mills A->data = (void*)b; 14230567c835SRichard Tran Mills 14240567c835SRichard Tran Mills ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 14250567c835SRichard Tran Mills 1426f4259b30SLisandro Dalcin b->w = NULL; 142755b25c41SPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRMPI(ierr); 14280567c835SRichard Tran Mills if (size == 1) { 14290567c835SRichard Tran Mills ierr = PetscObjectChangeTypeName((PetscObject)A,MATSEQKAIJ);CHKERRQ(ierr); 14300567c835SRichard Tran Mills A->ops->destroy = MatDestroy_SeqKAIJ; 1431bb6fb833SRichard Tran Mills A->ops->mult = MatMult_SeqKAIJ; 1432bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_SeqKAIJ; 1433bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_SeqKAIJ; 14340567c835SRichard Tran Mills A->ops->getrow = MatGetRow_SeqKAIJ; 14350567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_SeqKAIJ; 14360567c835SRichard Tran Mills A->ops->sor = MatSOR_SeqKAIJ; 14370567c835SRichard Tran Mills } else { 14380567c835SRichard Tran Mills ierr = PetscObjectChangeTypeName((PetscObject)A,MATMPIKAIJ);CHKERRQ(ierr); 14390567c835SRichard Tran Mills A->ops->destroy = MatDestroy_MPIKAIJ; 1440bb6fb833SRichard Tran Mills A->ops->mult = MatMult_MPIKAIJ; 1441bb6fb833SRichard Tran Mills A->ops->multadd = MatMultAdd_MPIKAIJ; 1442bb6fb833SRichard Tran Mills A->ops->invertblockdiagonal = MatInvertBlockDiagonal_MPIKAIJ; 14430567c835SRichard Tran Mills A->ops->getrow = MatGetRow_MPIKAIJ; 14440567c835SRichard Tran Mills A->ops->restorerow = MatRestoreRow_MPIKAIJ; 14450567c835SRichard Tran Mills ierr = PetscObjectComposeFunction((PetscObject)A,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPIKAIJ);CHKERRQ(ierr); 14460567c835SRichard Tran Mills } 1447*06d61a37SPierre Jolivet A->ops->setup = MatSetUp_KAIJ; 1448*06d61a37SPierre Jolivet A->ops->view = MatView_KAIJ; 14490567c835SRichard Tran Mills A->ops->createsubmatrix = MatCreateSubMatrix_KAIJ; 145049bd79ccSDebojyoti Ghosh PetscFunctionReturn(0); 145149bd79ccSDebojyoti Ghosh } 1452