xref: /petsc/src/mat/tests/ex20.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 
2 static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C,A;
9   PetscInt       i,j,m = 5,n = 4,Ii,J;
10   PetscMPIInt    rank,size;
11   PetscScalar    v;
12   char           mtype[256];
13 
14   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16 
17   /* This example does not work correctly for np > 2 */
18   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
19   PetscCheck(size <= 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Use np <= 2");
20 
21   /* Create the matrix for the five point stencil, YET AGAIN */
22   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
23   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
24   PetscCall(MatSetFromOptions(C));
25   PetscCall(MatSetUp(C));
26   for (i=0; i<m; i++) {
27     for (j=2*rank; j<2*rank+2; j++) {
28       v = -1.0;  Ii = j + n*i;
29       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
30       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
31       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
32       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
33       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
34     }
35   }
36 
37   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
38   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
39   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
40   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
41   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
42   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
43 
44   PetscCall(PetscStrcpy(mtype,MATSAME));
45   PetscCall(PetscOptionsGetString(NULL,NULL,"-conv_mat_type",mtype,sizeof(mtype),NULL));
46   PetscCall(MatConvert(C,mtype,MAT_INITIAL_MATRIX,&A));
47   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
48   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
49   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
50   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_IMPL));
51   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
52   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
53 
54   /* Free data structures */
55   PetscCall(MatDestroy(&A));
56   PetscCall(MatDestroy(&C));
57 
58   PetscCall(PetscFinalize());
59   return 0;
60 }
61 
62 /*TEST
63 
64    test:
65       args: -conv_mat_type seqaij
66 
67 TEST*/
68