1*c4762a1bSJed Brown static const char help[] = "Test MatGetLocalSubMatrix() with multiple levels of nesting.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc, char *argv[]) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown PetscErrorCode ierr; 8*c4762a1bSJed Brown IS is0a,is0b,is0,is1,isl0a,isl0b,isl0,isl1; 9*c4762a1bSJed Brown Mat A,Aexplicit; 10*c4762a1bSJed Brown PetscBool usenest; 11*c4762a1bSJed Brown PetscMPIInt rank,size; 12*c4762a1bSJed Brown PetscInt i,j; 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 15*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 16*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown { 19*c4762a1bSJed Brown PetscInt ix0a[1],ix0b[1],ix0[2],ix1[1]; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ix0a[0] = rank*2+0; 22*c4762a1bSJed Brown ix0b[0] = rank*2+1; 23*c4762a1bSJed Brown ix0[0] = rank*3+0; ix0[1] = rank*3+1; 24*c4762a1bSJed Brown ix1[0] = rank*3+2; 25*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,1,ix0a,PETSC_COPY_VALUES,&is0a);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,1,ix0b,PETSC_COPY_VALUES,&is0b);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,2,ix0,PETSC_COPY_VALUES,&is0);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,1,ix1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 29*c4762a1bSJed Brown } 30*c4762a1bSJed Brown { 31*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,6,0,1,&isl0);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,3,0,1,&isl0a);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,3,3,1,&isl0b);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,3,6,1,&isl1);CHKERRQ(ierr); 35*c4762a1bSJed Brown } 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown usenest = PETSC_FALSE; 38*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-nest",&usenest,NULL);CHKERRQ(ierr); 39*c4762a1bSJed Brown if (usenest) { 40*c4762a1bSJed Brown ISLocalToGlobalMapping l2g; 41*c4762a1bSJed Brown PetscInt l2gind[3]; 42*c4762a1bSJed Brown Mat B[9]; 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown l2gind[0] = (rank-1+size)%size; l2gind[1] = rank; l2gind[2] = (rank+1)%size; 45*c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,l2gind,PETSC_COPY_VALUES,&l2g);CHKERRQ(ierr); 46*c4762a1bSJed Brown for (i=0; i<9; i++) { 47*c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,1,1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&B[i]);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatSetUp(B[i]);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatSetLocalToGlobalMapping(B[i],l2g,l2g);CHKERRQ(ierr); 50*c4762a1bSJed Brown } 51*c4762a1bSJed Brown { 52*c4762a1bSJed Brown IS isx[2]; 53*c4762a1bSJed Brown Mat Bx00[4],Bx01[2],Bx10[2]; 54*c4762a1bSJed Brown Mat B00,B01,B10; 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown isx[0] = is0a; isx[1] = is0b; 57*c4762a1bSJed Brown Bx00[0] = B[0]; Bx00[1] = B[1]; Bx00[2] = B[3]; Bx00[3] = B[4]; 58*c4762a1bSJed Brown Bx01[0] = B[2]; Bx01[1] = B[5]; 59*c4762a1bSJed Brown Bx10[0] = B[6]; Bx10[1] = B[7]; 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,2,isx,2,isx,Bx00,&B00);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatSetUp(B00);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,2,isx,1,NULL,Bx01,&B01);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatSetUp(B01);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,1,NULL,2,isx,Bx10,&B10);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatSetUp(B10);CHKERRQ(ierr); 67*c4762a1bSJed Brown { 68*c4762a1bSJed Brown Mat By[4]; 69*c4762a1bSJed Brown IS isy[2]; 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown By[0] = B00; By[1] = B01; By[2] = B10; By[3] = B[8]; 72*c4762a1bSJed Brown isy[0] = is0; isy[1] = is1; 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,2,isy,2,isy,By,&A);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown ierr = MatDestroy(&B00);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatDestroy(&B01);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatDestroy(&B10);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown for (i=0; i<9; i++) {ierr = MatDestroy(&B[i]);CHKERRQ(ierr);} 82*c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr); 83*c4762a1bSJed Brown } else { 84*c4762a1bSJed Brown ISLocalToGlobalMapping l2g; 85*c4762a1bSJed Brown PetscInt l2gind[9]; 86*c4762a1bSJed Brown for (i=0; i<3; i++) for (j=0; j<3; j++) l2gind[3*i+j] = ((rank-1+j+size) % size)*3 + i; 87*c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,9,l2gind,PETSC_COPY_VALUES,&l2g);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatCreateAIJ(PETSC_COMM_WORLD,3,3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatSetLocalToGlobalMapping(A,l2g,l2g);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&l2g);CHKERRQ(ierr); 91*c4762a1bSJed Brown } 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown { 94*c4762a1bSJed Brown Mat A00,A11,A0a0a,A0a0b; 95*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A,isl0,isl0,&A00);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A,isl1,isl1,&A11);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A00,isl0a,isl0a,&A0a0a);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A00,isl0a,isl0b,&A0a0b);CHKERRQ(ierr); 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown ierr = MatSetValueLocal(A0a0a,0,0,100*rank+1,ADD_VALUES);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = MatSetValueLocal(A0a0a,0,1,100*rank+2,ADD_VALUES);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = MatSetValueLocal(A0a0a,2,2,100*rank+9,ADD_VALUES);CHKERRQ(ierr); 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown ierr = MatSetValueLocal(A0a0b,1,1,100*rank+50+5,ADD_VALUES);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown ierr = MatSetValueLocal(A11,0,0,1000*(rank+1)+1,ADD_VALUES);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = MatSetValueLocal(A11,1,2,1000*(rank+1)+6,ADD_VALUES);CHKERRQ(ierr); 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A00,isl0a,isl0a,&A0a0a);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A00,isl0a,isl0b,&A0a0b);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A,isl0,isl0,&A00);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A,isl1,isl1,&A11);CHKERRQ(ierr); 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown ierr = MatComputeOperator(A,MATAIJ,&Aexplicit);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = MatView(Aexplicit,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = MatDestroy(&Aexplicit);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = ISDestroy(&is0a);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = ISDestroy(&is0b);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = ISDestroy(&is0);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = ISDestroy(&isl0a);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = ISDestroy(&isl0b);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = ISDestroy(&isl0);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = ISDestroy(&isl1);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = PetscFinalize(); 131*c4762a1bSJed Brown return ierr; 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /*TEST 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown test: 138*c4762a1bSJed Brown nsize: 3 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown test: 141*c4762a1bSJed Brown suffix: nest 142*c4762a1bSJed Brown nsize: 3 143*c4762a1bSJed Brown args: -nest 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown TEST*/ 146