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; 22*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 23*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 24c4762a1bSJed Brown M = 10; N = 10; Z = 10; 25c4762a1bSJed Brown dof = 10; 26c4762a1bSJed Brown 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_3D",&Test_3D,NULL)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL)); 31c4762a1bSJed Brown /* Set up distributed array for fine grid */ 32c4762a1bSJed Brown if (!Test_3D) { 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,1,NULL,NULL,&coarsedm)); 34c4762a1bSJed Brown } else { 35*5f80ce2aSJacob 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)); 36c4762a1bSJed Brown } 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(coarsedm)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(coarsedm)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(coarsedm,PetscObjectComm((PetscObject)coarsedm),&finedm)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /*------------------------------------------------------------*/ 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(finedm,MATAIJ)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(finedm,&A)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* set val=one to A */ 48c4762a1bSJed Brown if (size == 1) { 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 50c4762a1bSJed Brown if (flg) { 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A,&array)); 52c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(A,&array)); 54c4762a1bSJed Brown } 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 56c4762a1bSJed Brown } else { 57c4762a1bSJed Brown Mat AA,AB; 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 60c4762a1bSJed Brown if (flg) { 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AA,&array)); 62c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AA,&array)); 64c4762a1bSJed Brown } 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 67c4762a1bSJed Brown if (flg) { 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AB,&array)); 69c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AB,&array)); 71c4762a1bSJed Brown } 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(coarsedm,finedm,&P,NULL)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 78c4762a1bSJed Brown /*------------------------------*/ 79c20d7725SJed Brown /* (1) Developer API */ 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(A,P,NULL,&C)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(C,MATPRODUCT_PtAP)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetAlgorithm(C,"allatonce")); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFill(C,PETSC_DEFAULT)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(C)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(C)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(C)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(C)); /* Test reuse of symbolic C */ 88c20d7725SJed Brown 89c2f6fc52SHong Zhang { /* Test MatProductView() */ 90c2f6fc52SHong Zhang PetscViewer viewer; 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(comm,NULL, &viewer)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductView(C,viewer)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 96c2f6fc52SHong Zhang } 97c2f6fc52SHong Zhang 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,P,C,10,&flg)); 992c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatProduct_PtAP"); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 101c20d7725SJed Brown 102c20d7725SJed Brown /* (2) User API */ 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C)); 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; 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 114c4762a1bSJed Brown 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,P,C,10,&flg)); 1162c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP"); 117c4762a1bSJed Brown 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&finedm)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&coarsedm)); 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 136c2f6fc52SHong Zhang args: -M 10 -N 10 -Z 10 137c2f6fc52SHong Zhang output_file: output/ex89_2.out 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: allatonce_merged 141c4762a1bSJed Brown nsize: 4 1423e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -mat_product_algorithm allatonce_merged 143c2f6fc52SHong Zhang output_file: output/ex89_3.out 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c2f6fc52SHong Zhang suffix: nonscalable_3D 147c4762a1bSJed Brown nsize: 4 1483e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm nonscalable 149c2f6fc52SHong Zhang output_file: output/ex89_4.out 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: allatonce_merged_3D 153c4762a1bSJed Brown nsize: 4 1543e662e0bSHong Zhang args: -M 10 -M 5 -M 10 -test_3D 1 -mat_product_algorithm allatonce_merged 155c2f6fc52SHong Zhang output_file: output/ex89_3.out 156c2f6fc52SHong Zhang 157c2f6fc52SHong Zhang test: 158c2f6fc52SHong Zhang suffix: nonscalable 159c2f6fc52SHong Zhang nsize: 4 1603e662e0bSHong Zhang args: -M 10 -N 10 -Z 10 -mat_product_algorithm nonscalable 161c2f6fc52SHong Zhang output_file: output/ex89_5.out 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164