xref: /petsc/src/mat/tests/ex21.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20c4762a1bSJed Brown   n    = 2*size;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
269566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(C,5,NULL));
28c4762a1bSJed Brown   for (i=0; i<m; i++) {
29c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
30c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
319566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
329566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
339566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
349566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
359566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
409566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
419566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
429566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
439566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
46c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
479566063dSJacob Faibussowitsch     PetscCall(MatGetRow(C,i,&nz,&idx,&values));
489566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i));
49c4762a1bSJed Brown     for (j=0; j<nz; j++) {
509566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g  ",idx[j],(double)PetscRealPart(values[j])));
51c4762a1bSJed Brown     }
529566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n"));
539566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(C,i,&nz,&idx,&values));
54c4762a1bSJed Brown   }
559566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A));
589566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
599566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
619566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
66b122ec5aSJacob Faibussowitsch   return 0;
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown /*TEST
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       args: -mat_type seqaij
73c4762a1bSJed Brown 
74c4762a1bSJed Brown TEST*/
75