xref: /petsc/src/mat/tests/ex162.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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   PetscErrorCode       ierr;
11   PetscMPIInt          size;
12 
13   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
14   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
15   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
16 
17   ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr);
18   ierr = MatSetSizes(A,4,4,4,4);CHKERRQ(ierr);
19   ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
20   ierr = MatSeqAIJSetPreallocation(A,4,NULL);CHKERRQ(ierr);
21 
22   row = 0; coli[0] = 1; coli[1] = 3; vali[0] = 1.0; vali[1] = 2.0;
23   ierr = MatSetValues(A,1,&row,2,coli,vali,ADD_VALUES);CHKERRQ(ierr);
24 
25   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;
26   ierr = MatSetValues(A,1,&row,4,coli,vali,ADD_VALUES);CHKERRQ(ierr);
27 
28   row = 2; coli[0] = 0; coli[1] = 3; vali[0] = 7.0; vali[1] = 8.0;
29   ierr = MatSetValues(A,1,&row,2,coli,vali,ADD_VALUES);CHKERRQ(ierr);
30 
31   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
32   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
33   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
34 
35   ierr = MatShift(A,0.0);CHKERRQ(ierr);
36   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
37 
38   ierr = MatDestroy(&A);CHKERRQ(ierr);
39   ierr = PetscFinalize();
40   return ierr;
41 }
42 
43 /*TEST
44 
45    test:
46 
47 TEST*/
48