1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with entire matrix, modified from ex59.c."; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **args) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat C,A,Adup; 9*c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,rstart,rend; 10*c4762a1bSJed Brown PetscMPIInt size,rank; 11*c4762a1bSJed Brown PetscErrorCode ierr; 12*c4762a1bSJed Brown PetscScalar v; 13*c4762a1bSJed Brown IS isrow; 14*c4762a1bSJed Brown PetscBool detect_bug = PETSC_FALSE; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 17*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-detect_bug",&detect_bug);CHKERRQ(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 ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown /* 28*c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 29*c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 30*c4762a1bSJed Brown */ 31*c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 32*c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 33*c4762a1bSJed Brown for (j=0; j<m*n; j++) { 34*c4762a1bSJed Brown v = i + j + 1; 35*c4762a1bSJed Brown ierr = MatSetValues(C,1,&i,1,&j,&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_COMMON);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 44*c4762a1bSJed Brown /* 45*c4762a1bSJed Brown Generate a new matrix consisting every row and column of the original matrix 46*c4762a1bSJed Brown */ 47*c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 50*c4762a1bSJed Brown if (detect_bug && !rank) { 51*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,1,rstart,1,&isrow);CHKERRQ(ierr); 52*c4762a1bSJed Brown } else { 53*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,rend-rstart,rstart,1,&isrow);CHKERRQ(ierr); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown ierr = MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */ 58*c4762a1bSJed Brown if (!rank ) { 59*c4762a1bSJed Brown i = 0; j = 0; v = 100; 60*c4762a1bSJed Brown ierr = MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown ierr = MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown /* Test MatDuplicate */ 71*c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&Adup);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = MatDestroy(&Adup);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = PetscFinalize(); 78*c4762a1bSJed Brown return ierr; 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown /*TEST 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown test: 85*c4762a1bSJed Brown nsize: 2 86*c4762a1bSJed Brown filter: grep -v "Mat Object" 87*c4762a1bSJed Brown requires: !complex 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown test: 90*c4762a1bSJed Brown suffix: 2 91*c4762a1bSJed Brown nsize: 3 92*c4762a1bSJed Brown args: -detect_bug 93*c4762a1bSJed Brown filter: grep -v "Mat Object" 94*c4762a1bSJed Brown requires: !complex 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown TEST*/ 97