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 DM coarsedm,finedm; 8c4762a1bSJed Brown PetscMPIInt size,rank; 9c4762a1bSJed Brown PetscInt M,N,Z,i,nrows; 10c4762a1bSJed Brown PetscScalar one = 1.0; 11c4762a1bSJed Brown PetscReal fill=2.0; 12c4762a1bSJed Brown Mat A,P,C; 13c20d7725SJed Brown PetscScalar *array,alpha; 14c4762a1bSJed Brown PetscBool Test_3D=PETSC_FALSE,flg; 15c4762a1bSJed Brown const PetscInt *ia,*ja; 16c4762a1bSJed Brown PetscInt dof; 17c4762a1bSJed Brown MPI_Comm comm; 18c4762a1bSJed Brown 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 20c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 225f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 23c4762a1bSJed Brown M = 10; N = 10; Z = 10; 24c4762a1bSJed Brown dof = 10; 25c4762a1bSJed Brown 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL)); 30c4762a1bSJed Brown /* Set up distributed array for fine grid */ 31c4762a1bSJed Brown if (!Test_3D) { 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm)); 33c4762a1bSJed Brown } else { 345f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 35c4762a1bSJed Brown } 365f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(coarsedm)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(coarsedm)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /*------------------------------------------------------------*/ 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(finedm,MATAIJ)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(finedm,&A)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* set val=one to A */ 47c4762a1bSJed Brown if (size == 1) { 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 49c4762a1bSJed Brown if (flg) { 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A,&array)); 51c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(A,&array)); 53c4762a1bSJed Brown } 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 55c4762a1bSJed Brown } else { 56c4762a1bSJed Brown Mat AA,AB; 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 59c4762a1bSJed Brown if (flg) { 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AA,&array)); 61c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AA,&array)); 63c4762a1bSJed Brown } 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 66c4762a1bSJed Brown if (flg) { 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AB,&array)); 68c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AB,&array)); 70c4762a1bSJed Brown } 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(coarsedm,finedm,&P,NULL)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 77c4762a1bSJed Brown /*------------------------------*/ 78c20d7725SJed Brown /* (1) Developer API */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(A,P,NULL,&C)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(C,MATPRODUCT_PtAP)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetAlgorithm(C,"allatonce")); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFill(C,PETSC_DEFAULT)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(C)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(C)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(C)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(C)); /* Test reuse of symbolic C */ 87c20d7725SJed Brown 88c2f6fc52SHong Zhang { /* Test MatProductView() */ 89c2f6fc52SHong Zhang PetscViewer viewer; 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(comm,NULL, &viewer)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductView(C,viewer)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 95c2f6fc52SHong Zhang } 96c2f6fc52SHong Zhang 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,P,C,10,&flg)); 9828b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatProduct_PtAP"); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 100c20d7725SJed Brown 101c20d7725SJed Brown /* (2) User API */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C)); 103c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 104c4762a1bSJed Brown alpha=1.0; 105c4762a1bSJed Brown for (i=0; i<1; i++) { 106c4762a1bSJed Brown alpha -= 0.1; 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 113c4762a1bSJed Brown 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,P,C,10,&flg)); 11528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP"); 116c4762a1bSJed Brown 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&finedm)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&coarsedm)); 122*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 123*b122ec5aSJacob Faibussowitsch return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*TEST 127c4762a1bSJed Brown 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown args: -M 10 -N 10 -Z 10 130c4762a1bSJed Brown output_file: output/ex89_1.out 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown suffix: allatonce 134c4762a1bSJed Brown nsize: 4 135c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10 136c2f6fc52SHong Zhang output_file: output/ex89_2.out 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: allatonce_merged 140c4762a1bSJed Brown nsize: 4 1413e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -mat_product_algorithm allatonce_merged 142c2f6fc52SHong Zhang output_file: output/ex89_3.out 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c2f6fc52SHong Zhang suffix: nonscalable_3D 146c4762a1bSJed Brown nsize: 4 1473e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm nonscalable 148c2f6fc52SHong Zhang output_file: output/ex89_4.out 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: allatonce_merged_3D 152c4762a1bSJed Brown nsize: 4 1533e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm allatonce_merged 154c2f6fc52SHong Zhang output_file: output/ex89_3.out 155c2f6fc52SHong Zhang 156c2f6fc52SHong Zhang test: 157c2f6fc52SHong Zhang suffix: nonscalable 158c2f6fc52SHong Zhang nsize: 4 1593e662e0bSHong Zhang args: -M 10 -N 10 -Z 10 -mat_product_algorithm nonscalable 160c2f6fc52SHong Zhang output_file: output/ex89_5.out 161c4762a1bSJed Brown 162c4762a1bSJed Brown TEST*/ 163