1*c4762a1bSJed Brown static char help[] = "Tests MatPermute() in parallel.\n\n"; 2*c4762a1bSJed Brown /* Results: 3*c4762a1bSJed Brown Sequential: 4*c4762a1bSJed Brown - seqaij: correct permutation 5*c4762a1bSJed Brown - seqbaij: permutation not supported for this MATTYPE 6*c4762a1bSJed Brown - seqsbaij: permutation not supported for this MATTYPE 7*c4762a1bSJed Brown Parallel: 8*c4762a1bSJed Brown - mpiaij: correct permutation 9*c4762a1bSJed Brown - mpibaij: correct permutation 10*c4762a1bSJed Brown - mpisbaij: permutation not supported for this MATTYPE 11*c4762a1bSJed Brown */ 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown #include <petscmat.h> 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown int main(int argc,char **argv) 16*c4762a1bSJed Brown { 17*c4762a1bSJed Brown const struct {PetscInt i,j; PetscScalar v;} entries[] = {{0,3,1.},{1,2,2.},{2,1,3.},{2,5,4.},{3,0,5.},{3,6,6.},{4,1,7.},{4,4,8.}}; 18*c4762a1bSJed Brown const PetscInt ixrow[5] = {4,2,1,0,3},ixcol[7] = {5,3,6,1,2,0,4}; 19*c4762a1bSJed Brown Mat A,B; 20*c4762a1bSJed Brown PetscErrorCode ierr; 21*c4762a1bSJed Brown PetscInt i,rstart,rend,cstart,cend; 22*c4762a1bSJed Brown IS isrow,iscol; 23*c4762a1bSJed Brown PetscViewer viewer; 24*c4762a1bSJed Brown PetscBool view_sparse; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 27*c4762a1bSJed Brown /* ------- Assemble matrix, --------- */ 28*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,7);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown for (i=0; i<(PetscInt)(sizeof(entries)/sizeof(entries[0])); i++) { 36*c4762a1bSJed Brown ierr = MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES);CHKERRQ(ierr); 37*c4762a1bSJed Brown } 38*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown /* ------ Prepare index sets ------ */ 42*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_WORLD,cend-cstart,ixcol+cstart,PETSC_USE_POINTER,&iscol);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = ISSetPermutation(isrow);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = ISSetPermutation(iscol);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown ierr = PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); 48*c4762a1bSJed Brown view_sparse = PETSC_FALSE; 49*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL);CHKERRQ(ierr); 50*c4762a1bSJed Brown if (!view_sparse) { 51*c4762a1bSJed Brown ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE);CHKERRQ(ierr); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Original matrix\n");CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatView(A,viewer);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown ierr = MatPermute(A,isrow,iscol,&B);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Permuted matrix\n");CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = MatView(B,viewer);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown if (!view_sparse) { 61*c4762a1bSJed Brown ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 62*c4762a1bSJed Brown } 63*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Row permutation\n");CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = ISView(isrow,viewer);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Column permutation\n");CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = ISView(iscol,viewer);CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown /* Free data structures */ 69*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = PetscFinalize(); 75*c4762a1bSJed Brown return ierr; 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown /*TEST 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown build: 83*c4762a1bSJed Brown requires: !complex 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown test: 86*c4762a1bSJed Brown args: -view_sparse 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown test: 89*c4762a1bSJed Brown suffix: 2 90*c4762a1bSJed Brown nsize: 2 91*c4762a1bSJed Brown args: -view_sparse 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown test: 94*c4762a1bSJed Brown suffix: 2b 95*c4762a1bSJed Brown nsize: 2 96*c4762a1bSJed Brown args: -mat_type baij -view_sparse 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown test: 99*c4762a1bSJed Brown suffix: 3 100*c4762a1bSJed Brown nsize: 3 101*c4762a1bSJed Brown args: -view_sparse 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown test: 104*c4762a1bSJed Brown suffix: 3b 105*c4762a1bSJed Brown nsize: 3 106*c4762a1bSJed Brown args: -mat_type baij -view_sparse 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown TEST*/ 109