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 PetscInt i,j,rowlens[] = {2,3,1},cols[] = {0,2,0,1,2,2}; 18c4762a1bSJed Brown PetscBool flg; 19c4762a1bSJed Brown 20*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-baij",&flg)); 22c4762a1bSJed Brown if (flg) { 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_WORLD,1,3,3,0,rowlens,&A)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetColumnIndices(A,cols)); 25c4762a1bSJed Brown } else { 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,3,3,0,rowlens,&A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetColumnIndices(A,cols)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown 30c4762a1bSJed Brown i = 0; j = 0; v = 1.0; 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 32c4762a1bSJed Brown i = 0; j = 2; v = 3.0; 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown i = 1; j = 0; v = 1.0; 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 37c4762a1bSJed Brown i = 1; j = 1; v = 2.0; 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 39c4762a1bSJed Brown i = 1; j = 2; v = 3.0; 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown i = 2; j = 2; v = 3.0; 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 44c4762a1bSJed Brown 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 50*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 51*b122ec5aSJacob Faibussowitsch return 0; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /*TEST 55c4762a1bSJed Brown 56c4762a1bSJed Brown test: 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown suffix: 2 60c4762a1bSJed Brown args: -baij 61c4762a1bSJed Brown 62c4762a1bSJed Brown TEST*/ 63