1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Tests sequential and parallel MatMatMatMult() and MatPtAP(). Modified from ex96.c \n\ 3c4762a1bSJed Brown -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\ 4c4762a1bSJed Brown -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\ 5c4762a1bSJed Brown -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\ 6c4762a1bSJed Brown -Npx <npx>, where <npx> = number of processors in the x-direction\n\ 7c4762a1bSJed Brown -Npy <npy>, where <npy> = number of processors in the y-direction\n\ 8c4762a1bSJed Brown -Npz <npz>, where <npz> = number of processors in the z-direction\n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Example of usage: mpiexec -n 3 ./ex41 -Mx 10 -My 10 -Mz 10 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscdm.h> 15c4762a1bSJed Brown #include <petscdmda.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* User-defined application contexts */ 18c4762a1bSJed Brown typedef struct { 19c4762a1bSJed Brown PetscInt mx,my,mz; /* number grid points in x, y and z direction */ 20c4762a1bSJed Brown Vec localX,localF; /* local vectors with ghost region */ 21c4762a1bSJed Brown DM da; 22c4762a1bSJed Brown Vec x,b,r; /* global vectors */ 23c4762a1bSJed Brown Mat J; /* Jacobian on grid */ 24c4762a1bSJed Brown } GridCtx; 25c4762a1bSJed Brown typedef struct { 26c4762a1bSJed Brown GridCtx fine; 27c4762a1bSJed Brown GridCtx coarse; 28c4762a1bSJed Brown PetscInt ratio; 29c4762a1bSJed Brown Mat Ii; /* interpolation from coarse to fine */ 30c4762a1bSJed Brown } AppCtx; 31c4762a1bSJed Brown 32c4762a1bSJed Brown #define COARSE_LEVEL 0 33c4762a1bSJed Brown #define FINE_LEVEL 1 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* 36c4762a1bSJed Brown Mm_ratio - ration of grid lines between fine and coarse grids. 37c4762a1bSJed Brown */ 38c4762a1bSJed Brown int main(int argc,char **argv) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown PetscErrorCode ierr; 41c4762a1bSJed Brown AppCtx user; 42c4762a1bSJed Brown PetscMPIInt size,rank; 43c4762a1bSJed Brown PetscInt m,n,M,N,i,nrows; 44c4762a1bSJed Brown PetscScalar one = 1.0; 45c4762a1bSJed Brown PetscReal fill=2.0; 46c20d7725SJed Brown Mat A,P,R,C,PtAP,D; 47c4762a1bSJed Brown PetscScalar *array; 48c4762a1bSJed Brown PetscRandom rdm; 49c4762a1bSJed Brown PetscBool Test_3D=PETSC_FALSE,flg; 50c4762a1bSJed Brown const PetscInt *ia,*ja; 51c4762a1bSJed Brown 52c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 535f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 545f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Get size of fine grids and coarse grids */ 57c4762a1bSJed Brown user.ratio = 2; 58c4762a1bSJed Brown user.coarse.mx = 4; user.coarse.my = 4; user.coarse.mz = 4; 59c4762a1bSJed Brown 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL)); 64c4762a1bSJed Brown if (user.coarse.mz) Test_3D = PETSC_TRUE; 65c4762a1bSJed Brown 66c4762a1bSJed Brown user.fine.mx = user.ratio*(user.coarse.mx-1)+1; 67c4762a1bSJed Brown user.fine.my = user.ratio*(user.coarse.my-1)+1; 68c4762a1bSJed Brown user.fine.mz = user.ratio*(user.coarse.mz-1)+1; 69c4762a1bSJed Brown 70dd400576SPatrick Sanan if (rank == 0) { 71c4762a1bSJed Brown if (!Test_3D) { 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"coarse grids: %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT "\n",user.coarse.mx,user.coarse.my,user.fine.mx,user.fine.my)); 73c4762a1bSJed Brown } else { 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"coarse grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "; fine grids: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",user.coarse.mx,user.coarse.my,user.coarse.mz,user.fine.mx,user.fine.my,user.fine.mz)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Set up distributed array for fine grid */ 79c4762a1bSJed Brown if (!Test_3D) { 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,user.fine.my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.fine.da)); 81c4762a1bSJed Brown } else { 82c4762a1bSJed Brown ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.fine.mx,user.fine.my,user.fine.mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, 83c4762a1bSJed Brown 1,1,NULL,NULL,NULL,&user.fine.da);CHKERRQ(ierr); 84c4762a1bSJed Brown } 855f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.fine.da)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.fine.da)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Create and set A at fine grids */ 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(user.fine.da,MATAIJ)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.fine.da,&A)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* set val=one to A (replace with random values!) */ 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 97c4762a1bSJed Brown if (size == 1) { 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 99c4762a1bSJed Brown if (flg) { 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A,&array)); 101c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(A,&array)); 103c4762a1bSJed Brown } 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 105c4762a1bSJed Brown } else { 106c4762a1bSJed Brown Mat AA,AB; 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 109c4762a1bSJed Brown if (flg) { 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AA,&array)); 111c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AA,&array)); 113c4762a1bSJed Brown } 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 116c4762a1bSJed Brown if (flg) { 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AB,&array)); 118c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AB,&array)); 120c4762a1bSJed Brown } 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown /* Set up distributed array for coarse grid */ 124c4762a1bSJed Brown if (!Test_3D) { 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,user.coarse.my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.coarse.da)); 126c4762a1bSJed Brown } else { 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,user.coarse.my,user.coarse.mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&user.coarse.da)); 128c4762a1bSJed Brown } 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.coarse.da)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.coarse.da)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 1335f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Get R = P^T */ 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(P,MAT_INITIAL_MATRIX,&R)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* C = R*A*P */ 139c20d7725SJed Brown /* Developer's API */ 1405f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(R,A,P,&D)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(D,MATPRODUCT_ABC)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(D)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(D)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(D)); /* Test reuse symbolic D */ 146c20d7725SJed Brown 147c20d7725SJed Brown /* User's API */ 148c20d7725SJed Brown { /* Test MatMatMatMult_Basic() */ 149c20d7725SJed Brown Mat Adense,Cdense; 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,Adense,P,MAT_INITIAL_MATRIX,fill,&Cdense)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,Adense,P,MAT_REUSE_MATRIX,fill,&Cdense)); 153c20d7725SJed Brown 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(D,Cdense,10,&flg)); 155*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D*v != Cdense*v"); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Adense)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 158c20d7725SJed Brown } 159c20d7725SJed Brown 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,fill,&C)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,fill,&C)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 163c4762a1bSJed Brown 164c20d7725SJed Brown /* Test D == C */ 1655f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(D,C,&flg)); 166*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D != C"); 167c20d7725SJed Brown 168c4762a1bSJed Brown /* Test C == PtAP */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&PtAP)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&PtAP)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(C,PtAP,&flg)); 172*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"C != PtAP"); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&PtAP)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* Clean up */ 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.fine.da)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.coarse.da)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 184c4762a1bSJed Brown ierr = PetscFinalize(); 185c4762a1bSJed Brown return ierr; 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /*TEST 189c4762a1bSJed Brown 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown 192c4762a1bSJed Brown test: 193c4762a1bSJed Brown suffix: 2 194c4762a1bSJed Brown nsize: 2 195c4762a1bSJed Brown args: -matmatmatmult_via scalable 196c4762a1bSJed Brown 197c4762a1bSJed Brown test: 198c4762a1bSJed Brown suffix: 3 199c4762a1bSJed Brown nsize: 2 200c4762a1bSJed Brown args: -matmatmatmult_via nonscalable 201c4762a1bSJed Brown output_file: output/ex111_1.out 202c4762a1bSJed Brown 203c4762a1bSJed Brown TEST*/ 204