1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\ 3*c4762a1bSJed Brown This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscmat.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc,char **args) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown Mat C,A; 10*c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,Ii,J,rstart,rend,nz; 11*c4762a1bSJed Brown PetscMPIInt rank,size; 12*c4762a1bSJed Brown PetscErrorCode ierr; 13*c4762a1bSJed Brown const PetscInt *idx; 14*c4762a1bSJed Brown PetscScalar v; 15*c4762a1bSJed Brown const PetscScalar *values; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 18*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 20*c4762a1bSJed Brown n = 2*size; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 23*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(C,5,NULL);CHKERRQ(ierr); 28*c4762a1bSJed Brown for (i=0; i<m; i++) { 29*c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 30*c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 31*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 32*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 33*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 34*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 35*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 36*c4762a1bSJed Brown } 37*c4762a1bSJed Brown } 38*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 46*c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 47*c4762a1bSJed Brown ierr = MatGetRow(C,i,&nz,&idx,&values);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %D: ",rank,i);CHKERRQ(ierr); 49*c4762a1bSJed Brown for (j=0; j<nz; j++) { 50*c4762a1bSJed Brown ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g ",idx[j],(double)PetscRealPart(values[j]));CHKERRQ(ierr); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown ierr = PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatRestoreRow(C,i,&nz,&idx,&values);CHKERRQ(ierr); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown ierr = MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscFinalize(); 66*c4762a1bSJed Brown return ierr; 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown /*TEST 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown test: 73*c4762a1bSJed Brown args: -mat_type seqaij 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown TEST*/ 76