xref: /petsc/src/mat/tests/ex21.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   n    = 2*size;
21 
22   /* create the matrix for the five point stencil, YET AGAIN*/
23   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
24   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
25   CHKERRQ(MatSetFromOptions(C));
26   CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
27   CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL));
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; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
32       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
33       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
34       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
35       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
36     }
37   }
38   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
39   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
40   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
41   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
42   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
43   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
44 
45   CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
46   for (i=rstart; i<rend; i++) {
47     CHKERRQ(MatGetRow(C,i,&nz,&idx,&values));
48     CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i));
49     for (j=0; j<nz; j++) {
50       CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g  ",idx[j],(double)PetscRealPart(values[j])));
51     }
52     CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n"));
53     CHKERRQ(MatRestoreRow(C,i,&nz,&idx,&values));
54   }
55   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
56 
57   CHKERRQ(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A));
58   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
59   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
60   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
61   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
62 
63   CHKERRQ(MatDestroy(&A));
64   CHKERRQ(MatDestroy(&C));
65   ierr = PetscFinalize();
66   return ierr;
67 }
68 
69 /*TEST
70 
71    test:
72       args: -mat_type seqaij
73 
74 TEST*/
75