xref: /petsc/src/mat/tests/ex20.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3) !
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   PetscFunctionBeginUser;
15   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
16   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
17 
18   /* This example does not work correctly for np > 2 */
19   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   PetscCheck(size <= 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Use np <= 2");
21 
22   /* Create the matrix for the five point stencil, YET AGAIN */
23   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
24   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
25   PetscCall(MatSetFromOptions(C));
26   PetscCall(MatSetUp(C));
27   for (i=0; i<m; i++) {
28     for (j=2*rank; j<2*rank+2; j++) {
29       v = -1.0;  Ii = j + n*i;
30       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
31       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
32       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
33       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
34       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
35     }
36   }
37 
38   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
39   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
40   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
41   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
42   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
43   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
44 
45   PetscCall(PetscStrcpy(mtype,MATSAME));
46   PetscCall(PetscOptionsGetString(NULL,NULL,"-conv_mat_type",mtype,sizeof(mtype),NULL));
47   PetscCall(MatConvert(C,mtype,MAT_INITIAL_MATRIX,&A));
48   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
49   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
50   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
51   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_IMPL));
52   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
53   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
54 
55   /* Free data structures */
56   PetscCall(MatDestroy(&A));
57   PetscCall(MatDestroy(&C));
58 
59   PetscCall(PetscFinalize());
60   return 0;
61 }
62 
63 /*TEST
64 
65    test:
66       args: -conv_mat_type seqaij
67 
68 TEST*/
69