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