1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] ="Tests sequential and parallel MatMatMatMult() and MatPtAP(). Modified from ex96.c \n\ 3*c4762a1bSJed Brown -Mx <xg>, where <xg> = number of coarse grid points in the x-direction\n\ 4*c4762a1bSJed Brown -My <yg>, where <yg> = number of coarse grid points in the y-direction\n\ 5*c4762a1bSJed Brown -Mz <zg>, where <zg> = number of coarse grid points in the z-direction\n\ 6*c4762a1bSJed Brown -Npx <npx>, where <npx> = number of processors in the x-direction\n\ 7*c4762a1bSJed Brown -Npy <npy>, where <npy> = number of processors in the y-direction\n\ 8*c4762a1bSJed Brown -Npz <npz>, where <npz> = number of processors in the z-direction\n\n"; 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown /* 11*c4762a1bSJed Brown Example of usage: mpiexec -n 3 ./ex41 -Mx 10 -My 10 -Mz 10 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown #include <petscdm.h> 15*c4762a1bSJed Brown #include <petscdmda.h> 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown /* User-defined application contexts */ 18*c4762a1bSJed Brown typedef struct { 19*c4762a1bSJed Brown PetscInt mx,my,mz; /* number grid points in x, y and z direction */ 20*c4762a1bSJed Brown Vec localX,localF; /* local vectors with ghost region */ 21*c4762a1bSJed Brown DM da; 22*c4762a1bSJed Brown Vec x,b,r; /* global vectors */ 23*c4762a1bSJed Brown Mat J; /* Jacobian on grid */ 24*c4762a1bSJed Brown } GridCtx; 25*c4762a1bSJed Brown typedef struct { 26*c4762a1bSJed Brown GridCtx fine; 27*c4762a1bSJed Brown GridCtx coarse; 28*c4762a1bSJed Brown PetscInt ratio; 29*c4762a1bSJed Brown Mat Ii; /* interpolation from coarse to fine */ 30*c4762a1bSJed Brown } AppCtx; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown #define COARSE_LEVEL 0 33*c4762a1bSJed Brown #define FINE_LEVEL 1 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* 36*c4762a1bSJed Brown Mm_ratio - ration of grid lines between fine and coarse grids. 37*c4762a1bSJed Brown */ 38*c4762a1bSJed Brown int main(int argc,char **argv) 39*c4762a1bSJed Brown { 40*c4762a1bSJed Brown PetscErrorCode ierr; 41*c4762a1bSJed Brown AppCtx user; 42*c4762a1bSJed Brown PetscMPIInt size,rank; 43*c4762a1bSJed Brown PetscInt m,n,M,N,i,nrows; 44*c4762a1bSJed Brown PetscScalar one = 1.0; 45*c4762a1bSJed Brown PetscReal fill=2.0; 46*c4762a1bSJed Brown Mat A,P,R,C,PtAP; 47*c4762a1bSJed Brown PetscScalar *array; 48*c4762a1bSJed Brown PetscRandom rdm; 49*c4762a1bSJed Brown PetscBool Test_3D=PETSC_FALSE,flg; 50*c4762a1bSJed Brown const PetscInt *ia,*ja; 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 53*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown /* Get size of fine grids and coarse grids */ 57*c4762a1bSJed Brown user.ratio = 2; 58*c4762a1bSJed Brown user.coarse.mx = 4; user.coarse.my = 4; user.coarse.mz = 4; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL);CHKERRQ(ierr); 64*c4762a1bSJed Brown if (user.coarse.mz) Test_3D = PETSC_TRUE; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown user.fine.mx = user.ratio*(user.coarse.mx-1)+1; 67*c4762a1bSJed Brown user.fine.my = user.ratio*(user.coarse.my-1)+1; 68*c4762a1bSJed Brown user.fine.mz = user.ratio*(user.coarse.mz-1)+1; 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown if (!rank) { 71*c4762a1bSJed Brown if (!Test_3D) { 72*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"coarse grids: %D %D; fine grids: %D %D\n",user.coarse.mx,user.coarse.my,user.fine.mx,user.fine.my);CHKERRQ(ierr); 73*c4762a1bSJed Brown } else { 74*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"coarse grids: %D %D %D; fine grids: %D %D %D\n",user.coarse.mx,user.coarse.my,user.coarse.mz,user.fine.mx,user.fine.my,user.fine.mz);CHKERRQ(ierr); 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown /* Set up distributed array for fine grid */ 79*c4762a1bSJed Brown if (!Test_3D) { 80*c4762a1bSJed 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); 81*c4762a1bSJed Brown } else { 82*c4762a1bSJed 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, 83*c4762a1bSJed Brown 1,1,NULL,NULL,NULL,&user.fine.da);CHKERRQ(ierr); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown ierr = DMSetFromOptions(user.fine.da);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = DMSetUp(user.fine.da);CHKERRQ(ierr); 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown /* Create and set A at fine grids */ 89*c4762a1bSJed Brown ierr = DMSetMatType(user.fine.da,MATAIJ);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = DMCreateMatrix(user.fine.da,&A);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown /* set val=one to A (replace with random values!) */ 95*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 97*c4762a1bSJed Brown if (size == 1) { 98*c4762a1bSJed Brown ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 99*c4762a1bSJed Brown if (flg) { 100*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr); 101*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 102*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr); 103*c4762a1bSJed Brown } 104*c4762a1bSJed Brown ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 105*c4762a1bSJed Brown } else { 106*c4762a1bSJed Brown Mat AA,AB; 107*c4762a1bSJed Brown ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 109*c4762a1bSJed Brown if (flg) { 110*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr); 111*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 112*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr); 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 116*c4762a1bSJed Brown if (flg) { 117*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr); 118*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 119*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr); 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown /* Set up distributed array for coarse grid */ 124*c4762a1bSJed Brown if (!Test_3D) { 125*c4762a1bSJed 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); 126*c4762a1bSJed Brown } else { 127*c4762a1bSJed 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); 128*c4762a1bSJed Brown } 129*c4762a1bSJed Brown ierr = DMSetFromOptions(user.coarse.da);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = DMSetUp(user.coarse.da);CHKERRQ(ierr); 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 133*c4762a1bSJed Brown ierr = DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL);CHKERRQ(ierr); 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /* Get R = P^T */ 136*c4762a1bSJed Brown ierr = MatTranspose(P,MAT_INITIAL_MATRIX,&R);CHKERRQ(ierr); 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown /* C = R*A*P */ 139*c4762a1bSJed Brown ierr = MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 142*c4762a1bSJed Brown 143*c4762a1bSJed Brown /* Test C == PtAP */ 144*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&PtAP);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&PtAP);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = MatEqual(C,PtAP,&flg);CHKERRQ(ierr); 147*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Matrices are not equal"); 148*c4762a1bSJed Brown ierr = MatDestroy(&PtAP);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* Clean up */ 151*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = DMDestroy(&user.fine.da);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = DMDestroy(&user.coarse.da);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatDestroy(&R);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = PetscFinalize(); 159*c4762a1bSJed Brown return ierr; 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown /*TEST 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown test: 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown test: 167*c4762a1bSJed Brown suffix: 2 168*c4762a1bSJed Brown nsize: 2 169*c4762a1bSJed Brown args: -matmatmatmult_via scalable 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown test: 172*c4762a1bSJed Brown suffix: 3 173*c4762a1bSJed Brown nsize: 2 174*c4762a1bSJed Brown args: -matmatmatmult_via nonscalable 175*c4762a1bSJed Brown output_file: output/ex111_1.out 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown TEST*/ 178