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