1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] ="Tests sequential and parallel DMCreateMatrix(), MatMatMult() and MatPtAP()\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 This test is modified from ~src/ksp/tests/ex19.c. 12c4762a1bSJed Brown Example of usage: mpiexec -n 3 ./ex96 -Mx 10 -My 10 -Mz 10 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* User-defined application contexts */ 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscInt mx,my,mz; /* number grid points in x, y and z direction */ 21c4762a1bSJed Brown Vec localX,localF; /* local vectors with ghost region */ 22c4762a1bSJed Brown DM da; 23c4762a1bSJed Brown Vec x,b,r; /* global vectors */ 24c4762a1bSJed Brown Mat J; /* Jacobian on grid */ 25c4762a1bSJed Brown } GridCtx; 26c4762a1bSJed Brown typedef struct { 27c4762a1bSJed Brown GridCtx fine; 28c4762a1bSJed Brown GridCtx coarse; 29c4762a1bSJed Brown PetscInt ratio; 30c4762a1bSJed Brown Mat Ii; /* interpolation from coarse to fine */ 31c4762a1bSJed Brown } AppCtx; 32c4762a1bSJed Brown 33c4762a1bSJed Brown #define COARSE_LEVEL 0 34c4762a1bSJed Brown #define FINE_LEVEL 1 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Mm_ratio - ration of grid lines between fine and coarse grids. 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown int main(int argc,char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown AppCtx user; 42c4762a1bSJed Brown PetscInt Npx=PETSC_DECIDE,Npy=PETSC_DECIDE,Npz=PETSC_DECIDE; 43c4762a1bSJed Brown PetscMPIInt size,rank; 44c4762a1bSJed Brown PetscInt m,n,M,N,i,nrows; 45c4762a1bSJed Brown PetscScalar one = 1.0; 46c4762a1bSJed Brown PetscReal fill=2.0; 47c4762a1bSJed Brown Mat A,A_tmp,P,C,C1,C2; 48c4762a1bSJed Brown PetscScalar *array,none = -1.0,alpha; 49c4762a1bSJed Brown Vec x,v1,v2,v3,v4; 50c4762a1bSJed Brown PetscReal norm,norm_tmp,norm_tmp1,tol=100.*PETSC_MACHINE_EPSILON; 51c4762a1bSJed Brown PetscRandom rdm; 52c4762a1bSJed Brown PetscBool Test_MatMatMult=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_3D=PETSC_TRUE,flg; 53c4762a1bSJed Brown const PetscInt *ia,*ja; 54c4762a1bSJed Brown 55*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tol",&tol,NULL)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown user.ratio = 2; 59c4762a1bSJed Brown user.coarse.mx = 20; user.coarse.my = 20; user.coarse.mz = 20; 60c4762a1bSJed Brown 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Mx",&user.coarse.mx,NULL)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-My",&user.coarse.my,NULL)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Mz",&user.coarse.mz,NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ratio",&user.ratio,NULL)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown if (user.coarse.mz) Test_3D = PETSC_TRUE; 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 695f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Npx",&Npx,NULL)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Npy",&Npy,NULL)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Npz",&Npz,NULL)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Set up distributed array for fine grid */ 75c4762a1bSJed Brown if (!Test_3D) { 765f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 77c4762a1bSJed Brown } else { 785f80ce2aSJacob 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,Npx,Npy,Npz,1,1,NULL,NULL,NULL,&user.coarse.da)); 79c4762a1bSJed Brown } 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.coarse.da)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.coarse.da)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* This makes sure the coarse DMDA has the same partition as the fine DMDA */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(user.coarse.da,PetscObjectComm((PetscObject)user.coarse.da),&user.fine.da)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Test DMCreateMatrix() */ 87c4762a1bSJed Brown /*------------------------------------------------------------*/ 885f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(user.fine.da,MATAIJ)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.fine.da,&A)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(user.fine.da,MATBAIJ)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.fine.da,&C)); 92c4762a1bSJed Brown 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATAIJ,MAT_INITIAL_MATRIX,&A_tmp)); /* not work for mpisbaij matrix! */ 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(A,A_tmp,&flg)); 9528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != C"); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A_tmp)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /*------------------------------------------------------------*/ 100c4762a1bSJed Brown 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 103dd400576SPatrick Sanan /* if (rank == 0) printf("A %d, %d\n",M,N); */ 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* set val=one to A */ 106c4762a1bSJed Brown if (size == 1) { 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 108c4762a1bSJed Brown if (flg) { 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(A,&array)); 110c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(A,&array)); 112c4762a1bSJed Brown } 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 114c4762a1bSJed Brown } else { 115c4762a1bSJed Brown Mat AA,AB; 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJGetSeqAIJ(A,&AA,&AB,NULL)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 118c4762a1bSJed Brown if (flg) { 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AA,&array)); 120c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AA,&array)); 122c4762a1bSJed Brown } 1235f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AA,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 125c4762a1bSJed Brown if (flg) { 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJGetArray(AB,&array)); 127c4762a1bSJed Brown for (i=0; i<ia[nrows]; i++) array[i] = one; 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJRestoreArray(AB,&array)); 129c4762a1bSJed Brown } 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(AB,0,PETSC_FALSE,PETSC_FALSE,&nrows,&ia,&ja,&flg)); 131c4762a1bSJed Brown } 1325f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); */ 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Create interpolation between the fine and coarse grids */ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(user.coarse.da,user.fine.da,&P,NULL)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(P,&m,&n)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(P,&M,&N)); 138dd400576SPatrick Sanan /* if (rank == 0) printf("P %d, %d\n",M,N); */ 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Create vectors v1 and v2 that are compatible with A */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&v1)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,NULL)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(v1,m,PETSC_DECIDE)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(v1)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v1,&v2)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Test MatMatMult(): C = A*P */ 150c4762a1bSJed Brown /*----------------------------*/ 151c4762a1bSJed Brown if (Test_MatMatMult) { 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A_tmp)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A_tmp,P,MAT_INITIAL_MATRIX,fill,&C)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 156c4762a1bSJed Brown alpha=1.0; 157c4762a1bSJed Brown for (i=0; i<2; i++) { 158c4762a1bSJed Brown alpha -= 0.1; 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A_tmp,alpha)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A_tmp,P,MAT_REUSE_MATRIX,fill,&C)); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 1635f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Test MatDuplicate() */ 166c4762a1bSJed Brown /*----------------------------*/ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C1,MAT_COPY_VALUES,&C2)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C2)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* Create vector x that is compatible with P */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(P,NULL,&n)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown norm = 0.0; 179c4762a1bSJed Brown for (i=0; i<10; i++) { 1805f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(P,x,v1)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A_tmp,v1,v2)); /* v2 = A*P*x */ 1835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,v1)); /* v1 = C*x */ 1845f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v1,none,v2)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v1,NORM_1,&norm_tmp)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v2,NORM_1,&norm_tmp1)); 187c4762a1bSJed Brown norm_tmp /= norm_tmp1; 188c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 189c4762a1bSJed Brown } 190dd400576SPatrick Sanan if (norm >= tol && rank == 0) { 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatMatMult(), |v1 - v2|/|v2|: %g\n",(double)norm)); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A_tmp)); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Test P^T * A * P - MatPtAP() */ 200c4762a1bSJed Brown /*------------------------------*/ 201c4762a1bSJed Brown if (Test_MatPtAP) { 2025f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(C,&m,&n)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 206c4762a1bSJed Brown alpha=1.0; 207c4762a1bSJed Brown for (i=0; i<1; i++) { 208c4762a1bSJed Brown alpha -= 0.1; 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C)); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* Free intermediate data structures created for reuse of C=Pt*A*P */ 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* Test MatDuplicate() */ 217c4762a1bSJed Brown /*----------------------------*/ 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C1,MAT_COPY_VALUES,&C2)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C2)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* Create vector x that is compatible with P */ 2245f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(P,&m,&n)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 228c4762a1bSJed Brown 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&v3)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(v3,n,PETSC_DECIDE)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(v3)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v3,&v4)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown norm = 0.0; 235c4762a1bSJed Brown for (i=0; i<10; i++) { 2365f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(P,x,v1)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v1,v2)); /* v2 = A*P*x */ 239c4762a1bSJed Brown 2405f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(P,v2,v3)); /* v3 = Pt*A*P*x */ 2415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,v4)); /* v3 = C*x */ 2425f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v4,none,v3)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v4,NORM_1,&norm_tmp)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v3,NORM_1,&norm_tmp1)); 245c4762a1bSJed Brown 246c4762a1bSJed Brown norm_tmp /= norm_tmp1; 247c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 248c4762a1bSJed Brown } 249dd400576SPatrick Sanan if (norm >= tol && rank == 0) { 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Error: MatPtAP(), |v3 - v4|/|v3|: %g\n",(double)norm)); 251c4762a1bSJed Brown } 2525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v3)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v4)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* Clean up */ 2595f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v1)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v2)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.fine.da)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.coarse.da)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 266*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 267*b122ec5aSJacob Faibussowitsch return 0; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown 270c4762a1bSJed Brown /*TEST 271c4762a1bSJed Brown 272c4762a1bSJed Brown test: 273c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 274c4762a1bSJed Brown output_file: output/ex96_1.out 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: nonscalable 278c4762a1bSJed Brown nsize: 3 279c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 280c4762a1bSJed Brown output_file: output/ex96_1.out 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: scalable 284c4762a1bSJed Brown nsize: 3 285c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable 286c4762a1bSJed Brown output_file: output/ex96_1.out 287c4762a1bSJed Brown 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown suffix: seq_scalable 290c4762a1bSJed Brown nsize: 3 2913e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm scalable -inner_offdiag_mat_product_algorithm scalable 292c4762a1bSJed Brown output_file: output/ex96_1.out 293c4762a1bSJed Brown 294c4762a1bSJed Brown test: 295c4762a1bSJed Brown suffix: seq_sorted 296c4762a1bSJed Brown nsize: 3 2973e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm sorted -inner_offdiag_mat_product_algorithm sorted 298c4762a1bSJed Brown output_file: output/ex96_1.out 299c4762a1bSJed Brown 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: seq_scalable_fast 302c4762a1bSJed Brown nsize: 3 3033e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm scalable_fast -inner_offdiag_mat_product_algorithm scalable_fast 304c4762a1bSJed Brown output_file: output/ex96_1.out 305c4762a1bSJed Brown 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown suffix: seq_heap 308c4762a1bSJed Brown nsize: 3 3093e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm heap -inner_offdiag_mat_product_algorithm heap 310c4762a1bSJed Brown output_file: output/ex96_1.out 311c4762a1bSJed Brown 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: seq_btheap 314c4762a1bSJed Brown nsize: 3 3153e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm btheap -inner_offdiag_mat_product_algorithm btheap 316c4762a1bSJed Brown output_file: output/ex96_1.out 317c4762a1bSJed Brown 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: seq_llcondensed 320c4762a1bSJed Brown nsize: 3 3213e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm llcondensed -inner_offdiag_mat_product_algorithm llcondensed 322c4762a1bSJed Brown output_file: output/ex96_1.out 323c4762a1bSJed Brown 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: seq_rowmerge 326c4762a1bSJed Brown nsize: 3 3273e662e0bSHong Zhang args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via scalable -inner_diag_mat_product_algorithm rowmerge -inner_offdiag_mat_product_algorithm rowmerge 328c4762a1bSJed Brown output_file: output/ex96_1.out 329c4762a1bSJed Brown 330c4762a1bSJed Brown test: 331c4762a1bSJed Brown suffix: allatonce 332c4762a1bSJed Brown nsize: 3 333c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce 334c4762a1bSJed Brown output_file: output/ex96_1.out 335c4762a1bSJed Brown 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: allatonce_merged 338c4762a1bSJed Brown nsize: 3 339c4762a1bSJed Brown args: -Mx 10 -My 5 -Mz 10 -matmatmult_via scalable -matptap_via allatonce_merged 340c4762a1bSJed Brown output_file: output/ex96_1.out 341c4762a1bSJed Brown 342c4762a1bSJed Brown TEST*/ 343