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