xref: /petsc/src/mat/tests/ex21.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\
3c4762a1bSJed Brown  This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat               C,A;
10c4762a1bSJed Brown   PetscInt          i,j,m = 3,n = 2,Ii,J,rstart,rend,nz;
11c4762a1bSJed Brown   PetscMPIInt       rank,size;
12c4762a1bSJed Brown   const PetscInt    *idx;
13c4762a1bSJed Brown   PetscScalar       v;
14c4762a1bSJed Brown   const PetscScalar *values;
15c4762a1bSJed Brown 
16*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
185f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
19c4762a1bSJed Brown   n    = 2*size;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL));
27c4762a1bSJed Brown   for (i=0; i<m; i++) {
28c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
29c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
305f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
315f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
325f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
335f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
345f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
35c4762a1bSJed Brown     }
36c4762a1bSJed Brown   }
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
43c4762a1bSJed Brown 
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend));
45c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
465f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetRow(C,i,&nz,&idx,&values));
475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i));
48c4762a1bSJed Brown     for (j=0; j<nz; j++) {
495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g  ",idx[j],(double)PetscRealPart(values[j])));
50c4762a1bSJed Brown     }
515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n"));
525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreRow(C,i,&nz,&idx,&values));
53c4762a1bSJed Brown   }
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
64*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
65*b122ec5aSJacob Faibussowitsch   return 0;
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*TEST
69c4762a1bSJed Brown 
70c4762a1bSJed Brown    test:
71c4762a1bSJed Brown       args: -mat_type seqaij
72c4762a1bSJed Brown 
73c4762a1bSJed Brown TEST*/
74