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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 17*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 18*9566063dSJacob Faibussowitsch PetscCallMPI(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*/ 22*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 23*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 24*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 25*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL)); 26*9566063dSJacob Faibussowitsch PetscCall(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; 30*9566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 31*9566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 32*9566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 33*9566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 34*9566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown } 37*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 38*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 39*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 40*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 41*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 42*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 43c4762a1bSJed Brown 44*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 45c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 46*9566063dSJacob Faibussowitsch PetscCall(MatGetRow(C,i,&nz,&idx,&values)); 47*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i)); 48c4762a1bSJed Brown for (j=0; j<nz; j++) { 49*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g ",idx[j],(double)PetscRealPart(values[j]))); 50c4762a1bSJed Brown } 51*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n")); 52*9566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(C,i,&nz,&idx,&values)); 53c4762a1bSJed Brown } 54*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout)); 55c4762a1bSJed Brown 56*9566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A)); 57*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 58*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 59*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 60*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 61c4762a1bSJed Brown 62*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 63*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 64*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 65b122ec5aSJacob 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