xref: /petsc/src/mat/tests/ex13.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            C,A;
9*c4762a1bSJed Brown   PetscInt       i,j,m = 5,n = 5,Ii,J;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscScalar    v;
12*c4762a1bSJed Brown   IS             perm,iperm;
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15*c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C);CHKERRQ(ierr);
16*c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
17*c4762a1bSJed Brown   for (i=0; i<m; i++) {
18*c4762a1bSJed Brown     for (j=0; j<n; j++) {
19*c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
20*c4762a1bSJed Brown       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
21*c4762a1bSJed Brown       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
22*c4762a1bSJed Brown       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
23*c4762a1bSJed Brown       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
24*c4762a1bSJed Brown       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
25*c4762a1bSJed Brown     }
26*c4762a1bSJed Brown   }
27*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown   ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = ISView(perm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = ISView(iperm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = ISDestroy(&iperm);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = PetscFinalize();
42*c4762a1bSJed Brown   return ierr;
43*c4762a1bSJed Brown }
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown /*TEST
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown    test:
49*c4762a1bSJed Brown       filter: grep -v "MPI processes"
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown TEST*/
52