xref: /petsc/src/mat/tests/ex21.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 
2 static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\
3  This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n";
4 
5 #include <petscmat.h>
6 
7 int main(int argc,char **args)
8 {
9   Mat               C,A;
10   PetscInt          i,j,m = 3,n = 2,Ii,J,rstart,rend,nz;
11   PetscMPIInt       rank,size;
12   PetscErrorCode    ierr;
13   const PetscInt    *idx;
14   PetscScalar       v;
15   const PetscScalar *values;
16 
17   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
19   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20   n    = 2*size;
21 
22   /* create the matrix for the five point stencil, YET AGAIN*/
23   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
24   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
25   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
26   ierr = MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);CHKERRQ(ierr);
27   ierr = MatSeqAIJSetPreallocation(C,5,NULL);CHKERRQ(ierr);
28   for (i=0; i<m; i++) {
29     for (j=2*rank; j<2*rank+2; j++) {
30       v = -1.0;  Ii = j + n*i;
31       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
32       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
33       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
34       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
35       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
36     }
37   }
38   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
39   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
41   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
44 
45   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
46   for (i=rstart; i<rend; i++) {
47     ierr = MatGetRow(C,i,&nz,&idx,&values);CHKERRQ(ierr);
48     ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %D: ",rank,i);CHKERRQ(ierr);
49     for (j=0; j<nz; j++) {
50       ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g  ",idx[j],(double)PetscRealPart(values[j]));CHKERRQ(ierr);
51     }
52     ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");CHKERRQ(ierr);
53     ierr = MatRestoreRow(C,i,&nz,&idx,&values);CHKERRQ(ierr);
54   }
55   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr);
56 
57   ierr = MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
58   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
59   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
60   ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
61   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
62 
63   ierr = MatDestroy(&A);CHKERRQ(ierr);
64   ierr = MatDestroy(&C);CHKERRQ(ierr);
65   ierr = PetscFinalize();
66   return ierr;
67 }
68 
69 
70 /*TEST
71 
72    test:
73       args: -mat_type seqaij
74 
75 TEST*/
76