xref: /petsc/src/mat/tests/ex13.c (revision 8cc725e69398de546bdc828d7b714aa2223f5218)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C,A;
9c4762a1bSJed Brown   PetscInt       i,j,m = 5,n = 5,Ii,J;
10c4762a1bSJed Brown   PetscScalar    v;
11c4762a1bSJed Brown   IS             perm,iperm;
12c4762a1bSJed Brown 
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
149566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C));
15c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
16c4762a1bSJed Brown   for (i=0; i<m; i++) {
17c4762a1bSJed Brown     for (j=0; j<n; j++) {
18c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
199566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
209566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
219566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
229566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
239566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
24c4762a1bSJed Brown     }
25c4762a1bSJed Brown   }
269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
28c4762a1bSJed Brown 
299566063dSJacob Faibussowitsch   PetscCall(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A));
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A,MATORDERINGND,&perm,&iperm));
329566063dSJacob Faibussowitsch   PetscCall(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
339566063dSJacob Faibussowitsch   PetscCall(ISView(iperm,PETSC_VIEWER_STDOUT_SELF));
349566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
409566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
41b122ec5aSJacob Faibussowitsch   return 0;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*TEST
45c4762a1bSJed Brown 
46c4762a1bSJed Brown    test:
47*8cc725e6SPierre Jolivet       filter: grep -v " MPI process"
48c4762a1bSJed Brown 
49c4762a1bSJed Brown TEST*/
50