1 static char help[] = "Tests MatShift for SeqAIJ matrices with some missing diagonal entries\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **argv) 6 { 7 Mat A; 8 PetscInt coli[4],row; 9 PetscScalar vali[4]; 10 PetscMPIInt size; 11 12 CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 13 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 14 PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 15 16 CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 17 CHKERRQ(MatSetSizes(A,4,4,4,4)); 18 CHKERRQ(MatSetType(A,MATSEQAIJ)); 19 CHKERRQ(MatSeqAIJSetPreallocation(A,4,NULL)); 20 21 row = 0; coli[0] = 1; coli[1] = 3; vali[0] = 1.0; vali[1] = 2.0; 22 CHKERRQ(MatSetValues(A,1,&row,2,coli,vali,ADD_VALUES)); 23 24 row = 1; coli[0] = 0; coli[1] = 1; coli[2] = 2; coli[3] = 3; vali[0] = 3.0; vali[1] = 4.0; vali[2] = 5.0; vali[3] = 6.0; 25 CHKERRQ(MatSetValues(A,1,&row,4,coli,vali,ADD_VALUES)); 26 27 row = 2; coli[0] = 0; coli[1] = 3; vali[0] = 7.0; vali[1] = 8.0; 28 CHKERRQ(MatSetValues(A,1,&row,2,coli,vali,ADD_VALUES)); 29 30 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 31 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 32 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 33 34 CHKERRQ(MatShift(A,0.0)); 35 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 36 37 CHKERRQ(MatDestroy(&A)); 38 CHKERRQ(PetscFinalize()); 39 return 0; 40 } 41 42 /*TEST 43 44 test: 45 46 TEST*/ 47