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; 53ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 54ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 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 60c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL);CHKERRQ(ierr); 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) { 728fffc762SJacob Faibussowitsch ierr = 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);CHKERRQ(ierr); 73c4762a1bSJed Brown } else { 748fffc762SJacob Faibussowitsch ierr = 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);CHKERRQ(ierr); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Set up distributed array for fine grid */ 79c4762a1bSJed Brown if (!Test_3D) { 80c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 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 } 85c4762a1bSJed Brown ierr = DMSetFromOptions(user.fine.da);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = DMSetUp(user.fine.da);CHKERRQ(ierr); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Create and set A at fine grids */ 89c4762a1bSJed Brown ierr = DMSetMatType(user.fine.da,MATAIJ);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = DMCreateMatrix(user.fine.da,&A);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* set val=one to A (replace with random values!) */ 95c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 97c4762a1bSJed Brown if (size == 1) { 98c4762a1bSJed Brown ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 99c4762a1bSJed Brown if (flg) { 100c4762a1bSJed Brown ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr); 101c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 102c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 105c4762a1bSJed Brown } else { 106c4762a1bSJed Brown Mat AA,AB; 107c4762a1bSJed Brown ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 109c4762a1bSJed Brown if (flg) { 110c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr); 111c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 112c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 116c4762a1bSJed Brown if (flg) { 117c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr); 118c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 119c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown /* Set up distributed array for coarse grid */ 124c4762a1bSJed Brown if (!Test_3D) { 125c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 126c4762a1bSJed Brown } else { 127c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown ierr = DMSetFromOptions(user.coarse.da);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = DMSetUp(user.coarse.da);CHKERRQ(ierr); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 133c4762a1bSJed Brown ierr = DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL);CHKERRQ(ierr); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Get R = P^T */ 136c4762a1bSJed Brown ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* C = R*A*P */ 139c20d7725SJed Brown /* Developer's API */ 140c20d7725SJed Brown ierr = MatProductCreate(R,A,P,&D);CHKERRQ(ierr); 141c20d7725SJed Brown ierr = MatProductSetType(D,MATPRODUCT_ABC);CHKERRQ(ierr); 142c20d7725SJed Brown ierr = MatProductSetFromOptions(D);CHKERRQ(ierr); 143c20d7725SJed Brown ierr = MatProductSymbolic(D);CHKERRQ(ierr); 144c20d7725SJed Brown ierr = MatProductNumeric(D);CHKERRQ(ierr); 145c20d7725SJed Brown ierr = MatProductNumeric(D);CHKERRQ(ierr); /* Test reuse symbolic D */ 146c20d7725SJed Brown 147c20d7725SJed Brown /* User's API */ 148c20d7725SJed Brown { /* Test MatMatMatMult_Basic() */ 149c20d7725SJed Brown Mat Adense,Cdense; 150c20d7725SJed Brown ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr); 151c20d7725SJed Brown ierr = MatMatMatMult(R,Adense,P,MAT_INITIAL_MATRIX,fill,&Cdense);CHKERRQ(ierr); 152c20d7725SJed Brown ierr = MatMatMatMult(R,Adense,P,MAT_REUSE_MATRIX,fill,&Cdense);CHKERRQ(ierr); 153c20d7725SJed Brown 154c20d7725SJed Brown ierr = MatMultEqual(D,Cdense,10,&flg);CHKERRQ(ierr); 155*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D*v != Cdense*v"); 156c20d7725SJed Brown ierr = MatDestroy(&Adense);CHKERRQ(ierr); 157c20d7725SJed Brown ierr = MatDestroy(&Cdense);CHKERRQ(ierr); 158c20d7725SJed Brown } 159c20d7725SJed Brown 160c4762a1bSJed Brown ierr = MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 1626718818eSStefano Zampini ierr = MatProductClear(C);CHKERRQ(ierr); 163c4762a1bSJed Brown 164c20d7725SJed Brown /* Test D == C */ 165c20d7725SJed Brown ierr = MatEqual(D,C,&flg);CHKERRQ(ierr); 166*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"D != C"); 167c20d7725SJed Brown 168c4762a1bSJed Brown /* Test C == PtAP */ 169c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&PtAP);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&PtAP);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatEqual(C,PtAP,&flg);CHKERRQ(ierr); 172*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"C != PtAP"); 173c4762a1bSJed Brown ierr = MatDestroy(&PtAP);CHKERRQ(ierr); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* Clean up */ 176c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = DMDestroy(&user.fine.da);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = DMDestroy(&user.coarse.da);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = MatDestroy(&R);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 183c20d7725SJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 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