1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatSeq(B)AIJSetColumnIndices().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /* 7c4762a1bSJed Brown Generate the following matrix: 8c4762a1bSJed Brown 9c4762a1bSJed Brown 1 0 3 10c4762a1bSJed Brown 1 2 3 11c4762a1bSJed Brown 0 0 3 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown int main(int argc,char **args) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown Mat A; 16c4762a1bSJed Brown PetscScalar v; 17c4762a1bSJed Brown PetscErrorCode ierr; 18c4762a1bSJed Brown PetscInt i,j,rowlens[] = {2,3,1},cols[] = {0,2,0,1,2,2}; 19c4762a1bSJed Brown PetscBool flg; 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-baij",&flg)); 23c4762a1bSJed Brown if (flg) { 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_WORLD,1,3,3,0,rowlens,&A)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetColumnIndices(A,cols)); 26c4762a1bSJed Brown } else { 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,3,3,0,rowlens,&A)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetColumnIndices(A,cols)); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31c4762a1bSJed Brown i = 0; j = 0; v = 1.0; 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 33c4762a1bSJed Brown i = 0; j = 2; v = 3.0; 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown i = 1; j = 0; v = 1.0; 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 38c4762a1bSJed Brown i = 1; j = 1; v = 2.0; 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 40c4762a1bSJed Brown i = 1; j = 2; v = 3.0; 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown i = 2; j = 2; v = 3.0; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 45c4762a1bSJed Brown 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 49c4762a1bSJed Brown 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 51c4762a1bSJed Brown ierr = PetscFinalize(); 52c4762a1bSJed Brown return ierr; 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /*TEST 56c4762a1bSJed Brown 57c4762a1bSJed Brown test: 58c4762a1bSJed Brown 59c4762a1bSJed Brown test: 60c4762a1bSJed Brown suffix: 2 61c4762a1bSJed Brown args: -baij 62c4762a1bSJed Brown 63c4762a1bSJed Brown TEST*/ 64