1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] ="Tests sequential and parallel DMCreateMatrix(), MatMatMult() and MatPtAP()\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 This test is modified from ~src/ksp/tests/ex19.c. 12*c4762a1bSJed Brown Example of usage: mpiexec -n 3 ./ex96 -Mx 10 -My 10 -Mz 10 13*c4762a1bSJed Brown */ 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown #include <petscdm.h> 16*c4762a1bSJed Brown #include <petscdmda.h> 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown /* User-defined application contexts */ 19*c4762a1bSJed Brown typedef struct { 20*c4762a1bSJed Brown PetscInt mx,my,mz; /* number grid points in x, y and z direction */ 21*c4762a1bSJed Brown Vec localX,localF; /* local vectors with ghost region */ 22*c4762a1bSJed Brown DM da; 23*c4762a1bSJed Brown Vec x,b,r; /* global vectors */ 24*c4762a1bSJed Brown Mat J; /* Jacobian on grid */ 25*c4762a1bSJed Brown } GridCtx; 26*c4762a1bSJed Brown typedef struct { 27*c4762a1bSJed Brown GridCtx fine; 28*c4762a1bSJed Brown GridCtx coarse; 29*c4762a1bSJed Brown PetscInt ratio; 30*c4762a1bSJed Brown Mat Ii; /* interpolation from coarse to fine */ 31*c4762a1bSJed Brown } AppCtx; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown #define COARSE_LEVEL 0 34*c4762a1bSJed Brown #define FINE_LEVEL 1 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown /* 37*c4762a1bSJed Brown Mm_ratio - ration of grid lines between fine and coarse grids. 38*c4762a1bSJed Brown */ 39*c4762a1bSJed Brown int main(int argc,char **argv) 40*c4762a1bSJed Brown { 41*c4762a1bSJed Brown PetscErrorCode ierr; 42*c4762a1bSJed Brown AppCtx user; 43*c4762a1bSJed Brown PetscInt Npx=PETSC_DECIDE,Npy=PETSC_DECIDE,Npz=PETSC_DECIDE; 44*c4762a1bSJed Brown PetscMPIInt size,rank; 45*c4762a1bSJed Brown PetscInt m,n,M,N,i,nrows; 46*c4762a1bSJed Brown PetscScalar one = 1.0; 47*c4762a1bSJed Brown PetscReal fill=2.0; 48*c4762a1bSJed Brown Mat A,A_tmp,P,C,C1,C2; 49*c4762a1bSJed Brown PetscScalar *array,none = -1.0,alpha; 50*c4762a1bSJed Brown Vec x,v1,v2,v3,v4; 51*c4762a1bSJed Brown PetscReal norm,norm_tmp,norm_tmp1,tol=100.*PETSC_MACHINE_EPSILON; 52*c4762a1bSJed Brown PetscRandom rdm; 53*c4762a1bSJed Brown PetscBool Test_MatMatMult=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_3D=PETSC_TRUE,flg; 54*c4762a1bSJed Brown const PetscInt *ia,*ja; 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 57*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown user.ratio = 2; 60*c4762a1bSJed Brown user.coarse.mx = 20; user.coarse.my = 20; user.coarse.mz = 20; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown if (user.coarse.mz) Test_3D = PETSC_TRUE; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Npx",&Npx,NULL);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Npy",&Npy,NULL);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-Npz",&Npz,NULL);CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown /* Set up distributed array for fine grid */ 76*c4762a1bSJed Brown if (!Test_3D) { 77*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,user.coarse.mx,user.coarse.my,Npx,Npy,1,1,NULL,NULL,&user.coarse.da);CHKERRQ(ierr); 78*c4762a1bSJed Brown } else { 79*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,Npx,Npy,Npz,1,1,NULL,NULL,NULL,&user.coarse.da);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown ierr = DMSetFromOptions(user.coarse.da);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = DMSetUp(user.coarse.da);CHKERRQ(ierr); 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 85*c4762a1bSJed Brown ierr = DMRefine(user.coarse.da,PetscObjectComm((PetscObject)user.coarse.da),&user.fine.da);CHKERRQ(ierr); 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown /* Test DMCreateMatrix() */ 88*c4762a1bSJed Brown /*------------------------------------------------------------*/ 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 = DMSetMatType(user.fine.da,MATBAIJ);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = DMCreateMatrix(user.fine.da,&C);CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown ierr = MatConvert(C,MATAIJ,MAT_INITIAL_MATRIX,&A_tmp);CHKERRQ(ierr); /* not work for mpisbaij matrix! */ 95*c4762a1bSJed Brown ierr = MatEqual(A,A_tmp,&flg);CHKERRQ(ierr); 96*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != C"); 97*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = MatDestroy(&A_tmp);CHKERRQ(ierr); 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown /*------------------------------------------------------------*/ 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 104*c4762a1bSJed Brown /* if (!rank) printf("A %d, %d\n",M,N); */ 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown /* set val=one to A */ 107*c4762a1bSJed Brown if (size == 1) { 108*c4762a1bSJed Brown ierr = MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 109*c4762a1bSJed Brown if (flg) { 110*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(A,&array);CHKERRQ(ierr); 111*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 112*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(A,&array);CHKERRQ(ierr); 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown ierr = MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 115*c4762a1bSJed Brown } else { 116*c4762a1bSJed Brown Mat AA,AB; 117*c4762a1bSJed Brown ierr = MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 119*c4762a1bSJed Brown if (flg) { 120*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AA,&array);CHKERRQ(ierr); 121*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 122*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AA,&array);CHKERRQ(ierr); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown ierr = MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 126*c4762a1bSJed Brown if (flg) { 127*c4762a1bSJed Brown ierr = MatSeqAIJGetArray(AB,&array);CHKERRQ(ierr); 128*c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 129*c4762a1bSJed Brown ierr = MatSeqAIJRestoreArray(AB,&array);CHKERRQ(ierr); 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown ierr = MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg);CHKERRQ(ierr); 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown /* ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 136*c4762a1bSJed Brown ierr = DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = MatGetLocalSize(P,&m,&n);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = MatGetSize(P,&M,&N);CHKERRQ(ierr); 139*c4762a1bSJed Brown /* if (!rank) printf("P %d, %d\n",M,N); */ 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* Create vectors v1 and v2 that are compatible with A */ 142*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&v1);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = VecSetSizes(v1,m,PETSC_DECIDE);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = VecSetFromOptions(v1);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = VecDuplicate(v1,&v2);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* Test MatMatMult(): C = A*P */ 151*c4762a1bSJed Brown /*----------------------------*/ 152*c4762a1bSJed Brown if (Test_MatMatMult) { 153*c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A_tmp);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = MatMatMult(A_tmp,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 157*c4762a1bSJed Brown alpha=1.0; 158*c4762a1bSJed Brown for (i=0; i<2; i++) { 159*c4762a1bSJed Brown alpha -= 0.1; 160*c4762a1bSJed Brown ierr = MatScale(A_tmp,alpha);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = MatMatMult(A_tmp,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 162*c4762a1bSJed Brown } 163*c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 164*c4762a1bSJed Brown ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 165*c4762a1bSJed Brown 166*c4762a1bSJed Brown /* Test MatDuplicate() */ 167*c4762a1bSJed Brown /*----------------------------*/ 168*c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_COPY_VALUES,&C1);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = MatDuplicate(C1,MAT_COPY_VALUES,&C2);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = MatDestroy(&C1);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = MatDestroy(&C2);CHKERRQ(ierr); 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown /* Create vector x that is compatible with P */ 174*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = MatGetLocalSize(P,NULL,&n);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown norm = 0.0; 180*c4762a1bSJed Brown for (i=0; i<10; i++) { 181*c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = MatMult(P,x,v1);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = MatMult(A_tmp,v1,v2);CHKERRQ(ierr); /* v2 = A*P*x */ 184*c4762a1bSJed Brown ierr = MatMult(C,x,v1);CHKERRQ(ierr); /* v1 = C*x */ 185*c4762a1bSJed Brown ierr = VecAXPY(v1,none,v2);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = VecNorm(v1,NORM_1,&norm_tmp);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = VecNorm(v2,NORM_1,&norm_tmp1);CHKERRQ(ierr); 188*c4762a1bSJed Brown norm_tmp /= norm_tmp1; 189*c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 190*c4762a1bSJed Brown } 191*c4762a1bSJed Brown if (norm >= tol && !rank) { 192*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult(), |v1 - v2|/|v2|: %g\n",(double)norm);CHKERRQ(ierr); 193*c4762a1bSJed Brown } 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = MatDestroy(&A_tmp);CHKERRQ(ierr); 198*c4762a1bSJed Brown } 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 201*c4762a1bSJed Brown /*------------------------------*/ 202*c4762a1bSJed Brown if (Test_MatPtAP) { 203*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = MatGetLocalSize(C,&m,&n);CHKERRQ(ierr); 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 207*c4762a1bSJed Brown alpha=1.0; 208*c4762a1bSJed Brown for (i=0; i<1; i++) { 209*c4762a1bSJed Brown alpha -= 0.1; 210*c4762a1bSJed Brown ierr = MatScale(A,alpha);CHKERRQ(ierr); 211*c4762a1bSJed Brown ierr = MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 215*c4762a1bSJed Brown ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 216*c4762a1bSJed Brown 217*c4762a1bSJed Brown /* Test MatDuplicate() */ 218*c4762a1bSJed Brown /*----------------------------*/ 219*c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_COPY_VALUES,&C1);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = MatDuplicate(C1,MAT_COPY_VALUES,&C2);CHKERRQ(ierr); 221*c4762a1bSJed Brown ierr = MatDestroy(&C1);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = MatDestroy(&C2);CHKERRQ(ierr); 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown /* Create vector x that is compatible with P */ 225*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = MatGetLocalSize(P,&m,&n);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&v3);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = VecSetSizes(v3,n,PETSC_DECIDE);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = VecSetFromOptions(v3);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = VecDuplicate(v3,&v4);CHKERRQ(ierr); 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown norm = 0.0; 236*c4762a1bSJed Brown for (i=0; i<10; i++) { 237*c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = MatMult(P,x,v1);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = MatMult(A,v1,v2);CHKERRQ(ierr); /* v2 = A*P*x */ 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown ierr = MatMultTranspose(P,v2,v3);CHKERRQ(ierr); /* v3 = Pt*A*P*x */ 242*c4762a1bSJed Brown ierr = MatMult(C,x,v4);CHKERRQ(ierr); /* v3 = C*x */ 243*c4762a1bSJed Brown ierr = VecAXPY(v4,none,v3);CHKERRQ(ierr); 244*c4762a1bSJed Brown ierr = VecNorm(v4,NORM_1,&norm_tmp);CHKERRQ(ierr); 245*c4762a1bSJed Brown ierr = VecNorm(v3,NORM_1,&norm_tmp1);CHKERRQ(ierr); 246*c4762a1bSJed Brown 247*c4762a1bSJed Brown norm_tmp /= norm_tmp1; 248*c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 249*c4762a1bSJed Brown } 250*c4762a1bSJed Brown if (norm >= tol && !rank) { 251*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",(double)norm);CHKERRQ(ierr); 252*c4762a1bSJed Brown } 253*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 254*c4762a1bSJed Brown ierr = VecDestroy(&v3);CHKERRQ(ierr); 255*c4762a1bSJed Brown ierr = VecDestroy(&v4);CHKERRQ(ierr); 256*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 257*c4762a1bSJed Brown } 258*c4762a1bSJed Brown 259*c4762a1bSJed Brown /* Clean up */ 260*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = VecDestroy(&v1);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = VecDestroy(&v2);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = DMDestroy(&user.fine.da);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = DMDestroy(&user.coarse.da);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = PetscFinalize(); 268*c4762a1bSJed Brown return ierr; 269*c4762a1bSJed Brown } 270*c4762a1bSJed Brown 271*c4762a1bSJed Brown /*TEST 272*c4762a1bSJed Brown 273*c4762a1bSJed Brown test: 274*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 275*c4762a1bSJed Brown output_file: output/ex96_1.out 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown test: 278*c4762a1bSJed Brown suffix: nonscalable 279*c4762a1bSJed Brown nsize: 3 280*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 281*c4762a1bSJed Brown output_file: output/ex96_1.out 282*c4762a1bSJed Brown 283*c4762a1bSJed Brown test: 284*c4762a1bSJed Brown suffix: scalable 285*c4762a1bSJed Brown nsize: 3 286*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable 287*c4762a1bSJed Brown output_file: output/ex96_1.out 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown test: 290*c4762a1bSJed Brown suffix: seq_scalable 291*c4762a1bSJed Brown nsize: 3 292*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via scalable -inner_offdiag_matmatmult_via scalable 293*c4762a1bSJed Brown output_file: output/ex96_1.out 294*c4762a1bSJed Brown 295*c4762a1bSJed Brown test: 296*c4762a1bSJed Brown suffix: seq_sorted 297*c4762a1bSJed Brown nsize: 3 298*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via sorted -inner_offdiag_matmatmult_via sorted 299*c4762a1bSJed Brown output_file: output/ex96_1.out 300*c4762a1bSJed Brown 301*c4762a1bSJed Brown test: 302*c4762a1bSJed Brown suffix: seq_scalable_fast 303*c4762a1bSJed Brown nsize: 3 304*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via scalable_fast -inner_offdiag_matmatmult_via scalable_fast 305*c4762a1bSJed Brown output_file: output/ex96_1.out 306*c4762a1bSJed Brown 307*c4762a1bSJed Brown test: 308*c4762a1bSJed Brown suffix: seq_heap 309*c4762a1bSJed Brown nsize: 3 310*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via heap -inner_offdiag_matmatmult_via heap 311*c4762a1bSJed Brown output_file: output/ex96_1.out 312*c4762a1bSJed Brown 313*c4762a1bSJed Brown test: 314*c4762a1bSJed Brown suffix: seq_btheap 315*c4762a1bSJed Brown nsize: 3 316*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via btheap -inner_offdiag_matmatmult_via btheap 317*c4762a1bSJed Brown output_file: output/ex96_1.out 318*c4762a1bSJed Brown 319*c4762a1bSJed Brown test: 320*c4762a1bSJed Brown suffix: seq_llcondensed 321*c4762a1bSJed Brown nsize: 3 322*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via llcondensed -inner_offdiag_matmatmult_via llcondensed 323*c4762a1bSJed Brown output_file: output/ex96_1.out 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown test: 326*c4762a1bSJed Brown suffix: seq_rowmerge 327*c4762a1bSJed Brown nsize: 3 328*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_matmatmult_via rowmerge -inner_offdiag_matmatmult_via rowmerge 329*c4762a1bSJed Brown output_file: output/ex96_1.out 330*c4762a1bSJed Brown 331*c4762a1bSJed Brown test: 332*c4762a1bSJed Brown suffix: allatonce 333*c4762a1bSJed Brown nsize: 3 334*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce 335*c4762a1bSJed Brown output_file: output/ex96_1.out 336*c4762a1bSJed Brown 337*c4762a1bSJed Brown test: 338*c4762a1bSJed Brown suffix: allatonce_merged 339*c4762a1bSJed Brown nsize: 3 340*c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce_merged 341*c4762a1bSJed Brown output_file: output/ex96_1.out 342*c4762a1bSJed Brown 343*c4762a1bSJed Brown TEST*/ 344