1c4762a1bSJed Brown static char help[] ="Tests MatPtAP() for MPIMAIJ and MPIAIJ \n "; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmda.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown PetscErrorCode ierr; 8c4762a1bSJed Brown DM coarsedm,finedm; 9c4762a1bSJed Brown PetscMPIInt size,rank; 10c4762a1bSJed Brown PetscInt M,N,Z,i,nrows; 11c4762a1bSJed Brown PetscScalar one = 1.0; 12c4762a1bSJed Brown PetscReal fill=2.0; 13c4762a1bSJed Brown Mat A,P,C; 14c20d7725SJed Brown PetscScalar *array,alpha; 15c4762a1bSJed Brown PetscBool Test_3D=PETSC_FALSE,flg; 16c4762a1bSJed Brown const PetscInt *ia,*ja; 17c4762a1bSJed Brown PetscInt dof; 18c4762a1bSJed Brown MPI_Comm comm; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 21c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 22ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 23ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 24c4762a1bSJed Brown M = 10; N = 10; Z = 10; 25c4762a1bSJed Brown dof = 10; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown /* Set up distributed array for fine grid */ 32c4762a1bSJed Brown if (!Test_3D) { 33c4762a1bSJed Brown ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm);CHKERRQ(ierr); 34c4762a1bSJed Brown } else { 35c4762a1bSJed Brown ierr = DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,Z,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,NULL,&coarsedm);CHKERRQ(ierr); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown ierr = DMSetFromOptions(coarsedm);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = DMSetUp(coarsedm);CHKERRQ(ierr); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 41c4762a1bSJed Brown ierr = DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm);CHKERRQ(ierr); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /*------------------------------------------------------------*/ 44c4762a1bSJed Brown ierr = DMSetMatType(finedm,MATAIJ);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = DMCreateMatrix(finedm,&A);CHKERRQ(ierr); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* set val=one to A */ 48c4762a1bSJed Brown if (size == 1) { 49c4762a1bSJed Brown ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 50c4762a1bSJed Brown if (flg) { 51c4762a1bSJed Brown ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr); 52c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 53c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 56c4762a1bSJed Brown } else { 57c4762a1bSJed Brown Mat AA,AB; 58c4762a1bSJed Brown ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 60c4762a1bSJed Brown if (flg) { 61c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr); 62c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 63c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 67c4762a1bSJed Brown if (flg) { 68c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr); 69c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 70c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 75c4762a1bSJed Brown ierr = DMCreateInterpolation(coarsedm,finedm,&P,NULL);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 78c4762a1bSJed Brown /*------------------------------*/ 79c20d7725SJed Brown /* (1) Developer API */ 80c20d7725SJed Brown ierr = MatProductCreate(A,P,NULL,&C);CHKERRQ(ierr); 81c20d7725SJed Brown ierr = MatProductSetType(C,MATPRODUCT_PtAP);CHKERRQ(ierr); 82*c2f6fc52SHong Zhang ierr = MatProductSetAlgorithm(C,"allatonce");CHKERRQ(ierr); 83c20d7725SJed Brown ierr = MatProductSetFill(C,PETSC_DEFAULT);CHKERRQ(ierr); 84c20d7725SJed Brown ierr = MatProductSetFromOptions(C);CHKERRQ(ierr); 85c20d7725SJed Brown ierr = MatProductSymbolic(C);CHKERRQ(ierr); 86c20d7725SJed Brown ierr = MatProductNumeric(C);CHKERRQ(ierr); 87c20d7725SJed Brown ierr = MatProductNumeric(C);CHKERRQ(ierr); /* Test reuse of symbolic C */ 88c20d7725SJed Brown 89*c2f6fc52SHong Zhang { /* Test MatProductView() */ 90*c2f6fc52SHong Zhang PetscViewer viewer; 91*c2f6fc52SHong Zhang ierr = PetscViewerASCIIOpen(comm,NULL, &viewer);CHKERRQ(ierr); 92*c2f6fc52SHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 93*c2f6fc52SHong Zhang ierr = MatProductView(C,viewer);CHKERRQ(ierr); 94*c2f6fc52SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 95*c2f6fc52SHong Zhang ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 96*c2f6fc52SHong Zhang } 97*c2f6fc52SHong Zhang 98c20d7725SJed Brown ierr = MatPtAPMultEqual(A,P,C,10,&flg);CHKERRQ(ierr); 99c20d7725SJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatProduct_PtAP"); 100c20d7725SJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 101c20d7725SJed Brown 102c20d7725SJed Brown /* (2) User API */ 103c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 104c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 105c4762a1bSJed Brown alpha=1.0; 106c4762a1bSJed Brown for (i=0; i<1; i++) { 107c4762a1bSJed Brown alpha -= 0.1; 108c4762a1bSJed Brown ierr = MatScale(A,alpha);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 1136718818eSStefano Zampini ierr = MatProductClear(C);CHKERRQ(ierr); 114c4762a1bSJed Brown 115c20d7725SJed Brown ierr = MatPtAPMultEqual(A,P,C,10,&flg);CHKERRQ(ierr); 116c20d7725SJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP"); 117c4762a1bSJed Brown 118c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMDestroy(&finedm);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = DMDestroy(&coarsedm);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = PetscFinalize(); 124c4762a1bSJed Brown return ierr; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*TEST 128c4762a1bSJed Brown 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown args: -M 10 -N 10 -Z 10 131c4762a1bSJed Brown output_file: output/ex89_1.out 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown suffix: allatonce 135c4762a1bSJed Brown nsize: 4 136*c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10 137*c2f6fc52SHong Zhang output_file: output/ex89_2.out 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: allatonce_merged 141c4762a1bSJed Brown nsize: 4 142*c2f6fc52SHong Zhang args: -M 10 -M 5 -M 10 -matproduct_ptap_via allatonce_merged 143*c2f6fc52SHong Zhang output_file: output/ex89_3.out 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146*c2f6fc52SHong Zhang suffix: nonscalable_3D 147c4762a1bSJed Brown nsize: 4 148*c2f6fc52SHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -matproduct_ptap_via nonscalable 149*c2f6fc52SHong Zhang output_file: output/ex89_4.out 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: allatonce_merged_3D 153c4762a1bSJed Brown nsize: 4 154*c2f6fc52SHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -matproduct_ptap_via allatonce_merged 155*c2f6fc52SHong Zhang output_file: output/ex89_3.out 156*c2f6fc52SHong Zhang 157*c2f6fc52SHong Zhang test: 158*c2f6fc52SHong Zhang suffix: nonscalable 159*c2f6fc52SHong Zhang nsize: 4 160*c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10 -matproduct_ptap_via nonscalable 161*c2f6fc52SHong Zhang output_file: output/ex89_5.out 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164