xref: /petsc/src/mat/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Creates a matrix, inserts some values, and tests MatCreateSubMatrices() and MatZeroEntries().\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **argv)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            mat,submat,submat1,*submatrices;
9*c4762a1bSJed Brown   PetscInt       m = 10,n = 10,i = 4,tmp,rstart,rend;
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   IS             irow,icol;
12*c4762a1bSJed Brown   PetscScalar    value = 1.0;
13*c4762a1bSJed Brown   PetscViewer    sviewer;
14*c4762a1bSJed Brown   PetscBool      allA = PETSC_FALSE;
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17*c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&mat);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = MatSetUp(mat);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
25*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
26*c4762a1bSJed Brown     value = (PetscReal)i+1; tmp = i % 5;
27*c4762a1bSJed Brown     ierr  = MatSetValues(mat,1,&tmp,1,&i,&value,INSERT_VALUES);CHKERRQ(ierr);
28*c4762a1bSJed Brown   }
29*c4762a1bSJed Brown   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original matrix\n");CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatView(mat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   /* Test MatCreateSubMatrix_XXX_All(), i.e., submatrix = A */
35*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-test_all",&allA,NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown   if (allA) {
37*c4762a1bSJed Brown     ierr   = ISCreateStride(PETSC_COMM_SELF,m,0,1,&irow);CHKERRQ(ierr);
38*c4762a1bSJed Brown     ierr   = ISCreateStride(PETSC_COMM_SELF,n,0,1,&icol);CHKERRQ(ierr);
39*c4762a1bSJed Brown     ierr   = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices);CHKERRQ(ierr);
40*c4762a1bSJed Brown     ierr   = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices);CHKERRQ(ierr);
41*c4762a1bSJed Brown     submat = *submatrices;
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown     /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */
44*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with all\n");CHKERRQ(ierr);
45*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"--------------------\n");CHKERRQ(ierr);
46*c4762a1bSJed Brown     ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
47*c4762a1bSJed Brown     ierr = MatView(submat,sviewer);CHKERRQ(ierr);
48*c4762a1bSJed Brown     ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
49*c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown     ierr = ISDestroy(&irow);CHKERRQ(ierr);
52*c4762a1bSJed Brown     ierr = ISDestroy(&icol);CHKERRQ(ierr);
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown     /* test getting a reference on a submat */
55*c4762a1bSJed Brown     ierr = PetscObjectReference((PetscObject)submat);CHKERRQ(ierr);
56*c4762a1bSJed Brown     ierr = MatDestroySubMatrices(1,&submatrices);CHKERRQ(ierr);
57*c4762a1bSJed Brown     ierr = MatDestroy(&submat);CHKERRQ(ierr);
58*c4762a1bSJed Brown   }
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /* Form submatrix with rows 2-4 and columns 4-8 */
61*c4762a1bSJed Brown   ierr   = ISCreateStride(PETSC_COMM_SELF,3,2,1,&irow);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr   = ISCreateStride(PETSC_COMM_SELF,5,4,1,&icol);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr   = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices);CHKERRQ(ierr);
64*c4762a1bSJed Brown   submat = *submatrices;
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   /* Test reuse submatrices */
67*c4762a1bSJed Brown   ierr = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices);CHKERRQ(ierr);
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown   /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */
70*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices\n");CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatView(submat,sviewer);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = PetscObjectReference((PetscObject)submat);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatDestroySubMatrices(1,&submatrices);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = MatDestroy(&submat);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* Form submatrix with rows 2-4 and all columns */
80*c4762a1bSJed Brown   ierr   = ISDestroy(&icol);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr   = ISCreateStride(PETSC_COMM_SELF,10,0,1,&icol);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr   = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr   = MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices);CHKERRQ(ierr);
84*c4762a1bSJed Brown   submat = *submatrices;
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with allcolumns\n");CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatView(submat,sviewer);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* Test MatDuplicate */
93*c4762a1bSJed Brown   ierr = MatDuplicate(submat,MAT_COPY_VALUES,&submat1);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = MatDestroy(&submat1);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Zero the original matrix */
97*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original zeroed matrix\n");CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatZeroEntries(mat);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = MatView(mat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   ierr = ISDestroy(&irow);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = ISDestroy(&icol);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscObjectReference((PetscObject)submat);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = MatDestroySubMatrices(1,&submatrices);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatDestroy(&submat);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatDestroy(&mat);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = PetscFinalize();
108*c4762a1bSJed Brown   return ierr;
109*c4762a1bSJed Brown }
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown /*TEST
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown    test:
116*c4762a1bSJed Brown       args: -mat_type aij
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown    test:
119*c4762a1bSJed Brown       suffix: 2
120*c4762a1bSJed Brown       args: -mat_type dense
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown    test:
123*c4762a1bSJed Brown       suffix: 3
124*c4762a1bSJed Brown       nsize: 3
125*c4762a1bSJed Brown       args: -mat_type aij
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown    test:
128*c4762a1bSJed Brown       suffix: 4
129*c4762a1bSJed Brown       nsize: 3
130*c4762a1bSJed Brown       args: -mat_type dense
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown    test:
133*c4762a1bSJed Brown       suffix: 5
134*c4762a1bSJed Brown       nsize: 3
135*c4762a1bSJed Brown       args: -mat_type aij -test_all
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown TEST*/
138