1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMult() and MatPtAP(), MatTransposeMatMult(), sequential MatMatTransposeMult(), MatRARt()\n\ 3c4762a1bSJed Brown Input arguments are:\n\ 4c4762a1bSJed Brown -f0 <input_file> -f1 <input_file> -f2 <input_file> -f3 <input_file> : file to load\n\n"; 5c4762a1bSJed Brown /* Example of usage: 6c4762a1bSJed Brown ./ex94 -f0 <A_binary> -f1 <B_binary> -matmatmult_mat_view ascii::ascii_info -matmatmulttr_mat_view 7c4762a1bSJed Brown mpiexec -n 3 ./ex94 -f0 medium -f1 medium -f2 arco1 -f3 arco1 -matmatmult_mat_view 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown B = A - B 14c4762a1bSJed Brown norm = norm(B) 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown PetscErrorCode MatNormDifference(Mat A,Mat B,PetscReal *norm) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscFunctionBegin; 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(B,NORM_FROBENIUS,norm)); 21c4762a1bSJed Brown PetscFunctionReturn(0); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24c4762a1bSJed Brown int main(int argc,char **args) 25c4762a1bSJed Brown { 26c4762a1bSJed Brown Mat A,A_save,B,AT,ATT,BT,BTT,P,R,C,C1; 27c4762a1bSJed Brown Vec x,v1,v2,v3,v4; 28c4762a1bSJed Brown PetscViewer viewer; 29c4762a1bSJed Brown PetscErrorCode ierr; 30c4762a1bSJed Brown PetscMPIInt size,rank; 31c4762a1bSJed Brown PetscInt i,m,n,j,*idxn,M,N,nzp,rstart,rend; 32c4762a1bSJed Brown PetscReal norm,norm_abs,norm_tmp,fill=4.0; 33c4762a1bSJed Brown PetscRandom rdm; 34c4762a1bSJed Brown char file[4][128]; 35c4762a1bSJed Brown PetscBool flg,preload = PETSC_TRUE; 36c4762a1bSJed Brown PetscScalar *a,rval,alpha,none = -1.0; 37c4762a1bSJed Brown PetscBool Test_MatMatMult=PETSC_TRUE,Test_MatMatTr=PETSC_TRUE,Test_MatPtAP=PETSC_TRUE,Test_MatRARt=PETSC_TRUE,Test_MatMatMatMult=PETSC_TRUE; 38c4762a1bSJed Brown PetscBool Test_MatAXPY=PETSC_FALSE,view=PETSC_FALSE; 39c4762a1bSJed Brown PetscInt pm,pn,pM,pN; 40c4762a1bSJed Brown MatInfo info; 41c4762a1bSJed Brown PetscBool seqaij; 42c4762a1bSJed Brown MatType mattype; 43c4762a1bSJed Brown Mat Cdensetest,Pdense,Cdense,Adense; 44c4762a1bSJed Brown 45c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 46*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 47*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-matops_view",&view,NULL)); 51c4762a1bSJed Brown if (view) { 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Load the matrices A_save and B */ 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f0",file[0],sizeof(file[0]),&flg)); 572c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for small matrix A with the -f0 option."); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f1",file[1],sizeof(file[1]),&flg)); 592c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for small matrix B with the -f1 option."); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f2",file[2],sizeof(file[2]),&flg)); 61c4762a1bSJed Brown if (!flg) { 62c4762a1bSJed Brown preload = PETSC_FALSE; 63c4762a1bSJed Brown } else { 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f3",file[3],sizeof(file[3]),&flg)); 652c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for test matrix B with the -f3 option."); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscPreLoadBegin(preload,"Load system"); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2*PetscPreLoadIt],FILE_MODE_READ,&viewer)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A_save)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A_save)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A_save,viewer)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 74c4762a1bSJed Brown 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2*PetscPreLoadIt+1],FILE_MODE_READ,&viewer)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(B,viewer)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 80c4762a1bSJed Brown 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(B,&mattype)); 82c4762a1bSJed Brown 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(B,&M,&N)); 84c4762a1bSJed Brown nzp = PetscMax((PetscInt)(0.1*M),5); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc((nzp+1)*(sizeof(PetscInt)+sizeof(PetscScalar)),&idxn)); 86c4762a1bSJed Brown a = (PetscScalar*)(idxn + nzp); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Create vectors v1 and v2 that are compatible with A_save */ 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&v1)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A_save,&m,NULL)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(v1,m,PETSC_DECIDE)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(v1)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v1,&v2)); 94c4762a1bSJed Brown 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* Test MatAXPY() */ 100c4762a1bSJed Brown /*-------------------*/ 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_MatAXPY",&Test_MatAXPY)); 102c4762a1bSJed Brown if (Test_MatAXPY) { 103c4762a1bSJed Brown Mat Btmp; 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A_save,MAT_COPY_VALUES,&A)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(B,MAT_COPY_VALUES,&Btmp)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(A,-1.0,B,DIFFERENT_NONZERO_PATTERN)); /* A = -B + A_save */ 107c4762a1bSJed Brown 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,-1.0)); /* A = -A = B - A_save */ 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Btmp,-1.0,A,DIFFERENT_NONZERO_PATTERN)); /* Btmp = -A + B = A_save */ 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A_save,Btmp,10,&flg)); 1112c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatAXPY() is incorrect"); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Btmp)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown Test_MatMatMult = PETSC_FALSE; 116c4762a1bSJed Brown Test_MatMatTr = PETSC_FALSE; 117c4762a1bSJed Brown Test_MatPtAP = PETSC_FALSE; 118c4762a1bSJed Brown Test_MatRARt = PETSC_FALSE; 119c4762a1bSJed Brown Test_MatMatMatMult = PETSC_FALSE; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* 1) Test MatMatMult() */ 123c4762a1bSJed Brown /* ---------------------*/ 124c4762a1bSJed Brown if (Test_MatMatMult) { 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A_save,MAT_COPY_VALUES,&A)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(A,&AT)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(AT,&ATT)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(B,&BT)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateTranspose(BT,&BTT)); 130c20d7725SJed Brown 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&C)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(AT,B,C,10,&flg)); 1332c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=AT*B"); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 135c20d7725SJed Brown 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(ATT,B,MAT_INITIAL_MATRIX,fill,&C)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(ATT,B,C,10,&flg)); 1382c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=ATT*B"); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 140c20d7725SJed Brown 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,C,10,&flg)); 1432c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for reuse C=A*B"); 144c20d7725SJed Brown /* ATT has different matrix type as A (although they have same internal data structure), 145c20d7725SJed Brown we cannot call MatProductReplaceMats(ATT,NULL,NULL,C) and MatMatMult(ATT,B,MAT_REUSE_MATRIX,fill,&C) */ 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 147c20d7725SJed Brown 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,BTT,MAT_INITIAL_MATRIX,fill,&C)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,BTT,C,10,&flg)); 1502c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult() for C=A*BTT"); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 152c20d7725SJed Brown 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(ATT,BTT,MAT_INITIAL_MATRIX,fill,&C)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,C,10,&flg)); 1552c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 157c20d7725SJed Brown 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&BTT)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&BT)); 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ATT)); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AT)); 162c20d7725SJed Brown 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C,"matmatmult_")); /* enable option '-matmatmult_' for matrix C */ 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 168c4762a1bSJed Brown alpha=1.0; 169c4762a1bSJed Brown for (i=0; i<2; i++) { 170c4762a1bSJed Brown alpha -=0.1; 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C)); 173c4762a1bSJed Brown } 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(A,B,C,10,&flg)); 1752c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatMatMult()"); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Test MatDuplicate() of C=A*B */ 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 182c4762a1bSJed Brown } /* if (Test_MatMatMult) */ 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* 2) Test MatTransposeMatMult() and MatMatTransposeMult() */ 185c4762a1bSJed Brown /* ------------------------------------------------------- */ 186c4762a1bSJed Brown if (Test_MatMatTr) { 187c4762a1bSJed Brown /* Create P */ 188c4762a1bSJed Brown PetscInt PN,rstart,rend; 189c4762a1bSJed Brown PN = M/2; 190c4762a1bSJed Brown nzp = 5; /* num of nonzeros in each row of P */ 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&P)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,M,PN)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(P,mattype)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(P,nzp,NULL)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(P,nzp,NULL,nzp,NULL)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(P,&rstart,&rend)); 197c4762a1bSJed Brown for (i=0; i<nzp; i++) { 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&a[i])); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 201c4762a1bSJed Brown for (j=0; j<nzp; j++) { 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 203c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval)*PN); 204c4762a1bSJed Brown } 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(P,1,&i,nzp,idxn,a,ADD_VALUES)); 206c4762a1bSJed Brown } 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* Create R = P^T */ 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(P,MAT_INITIAL_MATRIX,&R)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown { /* Test R = P^T, C1 = R*B */ 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(R,B,MAT_INITIAL_MATRIX,fill,&C1)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(P,MAT_REUSE_MATRIX,&R)); 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(R,B,MAT_REUSE_MATRIX,fill,&C1)); 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* C = P^T*B */ 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(P,B,MAT_INITIAL_MATRIX,fill,&C)); 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 225*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(P,B,MAT_REUSE_MATRIX,fill,&C)); 226c4762a1bSJed Brown if (view) { 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C = P^T * B:\n")); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 229c4762a1bSJed Brown } 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 231c4762a1bSJed Brown if (view) { 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * B after MatProductClear():\n")); 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Compare P^T*B and R*B */ 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(R,B,MAT_INITIAL_MATRIX,fill,&C1)); 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNormDifference(C,C1,&norm)); 2392c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatTransposeMatMult(): %g",(double)norm); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* Test MatDuplicate() of C=P^T*B */ 243*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* C = B*R^T */ 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij)); 249c4762a1bSJed Brown if (size == 1 && seqaij) { 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(B,R,MAT_INITIAL_MATRIX,fill,&C)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C,"matmatmulttr_")); /* enable '-matmatmulttr_' for matrix C */ 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(B,R,MAT_REUSE_MATRIX,fill,&C)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* Check */ 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,P,MAT_INITIAL_MATRIX,fill,&C1)); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNormDifference(C,C1,&norm)); 2602c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatMatTransposeMult() %g",(double)norm); 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C1)); 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 263c4762a1bSJed Brown } 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* 3) Test MatPtAP() */ 269c4762a1bSJed Brown /*-------------------*/ 270c4762a1bSJed Brown if (Test_MatPtAP) { 271c4762a1bSJed Brown PetscInt PN; 272c4762a1bSJed Brown Mat Cdup; 273c4762a1bSJed Brown 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A_save,MAT_COPY_VALUES,&A)); 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 277c4762a1bSJed Brown 278c4762a1bSJed Brown PN = M/2; 279c4762a1bSJed Brown nzp = (PetscInt)(0.1*PN+1); /* num of nozeros in each row of P */ 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&P)); 281*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,N,PN)); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(P,mattype)); 283*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(P,nzp,NULL)); 284*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(P,nzp,NULL,nzp,NULL)); 285c4762a1bSJed Brown for (i=0; i<nzp; i++) { 286*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&a[i])); 287c4762a1bSJed Brown } 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(P,&rstart,&rend)); 289c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 290c4762a1bSJed Brown for (j=0; j<nzp; j++) { 291*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 292c4762a1bSJed Brown idxn[j] = (PetscInt)(PetscRealPart(rval)*PN); 293c4762a1bSJed Brown } 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(P,1,&i,nzp,idxn,a,ADD_VALUES)); 295c4762a1bSJed Brown } 296*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 298c4762a1bSJed Brown 299*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatView(P,PETSC_VIEWER_STDOUT_WORLD)); */ 300*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(P,&pM,&pN)); 301*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(P,&pm,&pn)); 302*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_INITIAL_MATRIX,fill,&C)); 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX - reuse symbolic C */ 305c4762a1bSJed Brown alpha=1.0; 306c4762a1bSJed Brown for (i=0; i<2; i++) { 307c4762a1bSJed Brown alpha -=0.1; 308*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,alpha)); 309*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,P,MAT_REUSE_MATRIX,fill,&C)); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* Test PtAP ops with P Dense and A either AIJ or SeqDense (it assumes MatPtAP_XAIJ_XAIJ is fine) */ 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij)); 314c4762a1bSJed Brown if (seqaij) { 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATSEQDENSE,MAT_INITIAL_MATRIX,&Cdensetest)); 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(P,MATSEQDENSE,MAT_INITIAL_MATRIX,&Pdense)); 317c4762a1bSJed Brown } else { 318*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdensetest)); 319*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(P,MATMPIDENSE,MAT_INITIAL_MATRIX,&Pdense)); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c20d7725SJed Brown /* test with A(AIJ), Pdense -- call MatPtAP_Basic() when np>1 */ 323*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,Pdense,MAT_INITIAL_MATRIX,fill,&Cdense)); 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,Pdense,MAT_REUSE_MATRIX,fill,&Cdense)); 325*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,Pdense,Cdense,10,&flg)); 3262c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in MatPtAP with A AIJ and P Dense"); 327*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* test with A SeqDense */ 330c4762a1bSJed Brown if (seqaij) { 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense)); 332*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(Adense,Pdense,MAT_INITIAL_MATRIX,fill,&Cdense)); 333*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(Adense,Pdense,MAT_REUSE_MATRIX,fill,&Cdense)); 334*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(Adense,Pdense,Cdense,10,&flg)); 3352c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in MatPtAP with A SeqDense and P SeqDense"); 336*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Adense)); 338c4762a1bSJed Brown } 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdensetest)); 340*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Pdense)); 341c4762a1bSJed Brown 342c4762a1bSJed Brown /* Test MatDuplicate() of C=PtAP and MatView(Cdup,...) */ 343*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&Cdup)); 344c4762a1bSJed Brown if (view) { 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * A * P:\n")); 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 347c4762a1bSJed Brown 348*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductClear(C)); 349*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nC = P^T * A * P after MatProductClear():\n")); 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 351c4762a1bSJed Brown 352*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nCdup:\n")); 353*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Cdup,PETSC_VIEWER_STDOUT_WORLD)); 354c4762a1bSJed Brown } 355*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdup)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown if (size>1 || !seqaij) Test_MatRARt = PETSC_FALSE; 358c4762a1bSJed Brown /* 4) Test MatRARt() */ 359c4762a1bSJed Brown /* ----------------- */ 360c4762a1bSJed Brown if (Test_MatRARt) { 361c20d7725SJed Brown Mat R, RARt, Rdense, RARtdense; 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(P,MAT_INITIAL_MATRIX,&R)); 363c20d7725SJed Brown 364c20d7725SJed Brown /* Test MatRARt_Basic(), MatMatMatMult_Basic() */ 365*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(R,MATDENSE,MAT_INITIAL_MATRIX,&Rdense)); 366*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,Rdense,MAT_INITIAL_MATRIX,2.0,&RARtdense)); 367*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,Rdense,MAT_REUSE_MATRIX,2.0,&RARtdense)); 368c20d7725SJed Brown 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(RARtdense,MATAIJ,MAT_INITIAL_MATRIX,&RARt)); 370*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNormDifference(C,RARt,&norm)); 3712c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"|PtAP - RARtdense| = %g",(double)norm); 372*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Rdense)); 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RARtdense)); 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RARt)); 375c20d7725SJed Brown 376c20d7725SJed Brown /* Test MatRARt() for aij matrices */ 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,R,MAT_INITIAL_MATRIX,2.0,&RARt)); 378*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRARt(A,R,MAT_REUSE_MATRIX,2.0,&RARt)); 379*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNormDifference(C,RARt,&norm)); 3802c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"|PtAP - RARt| = %g",(double)norm); 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RARt)); 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (Test_MatMatMatMult && size == 1) { 386c4762a1bSJed Brown Mat R, RAP; 387*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(P,MAT_INITIAL_MATRIX,&R)); 388*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,A,P,MAT_INITIAL_MATRIX,2.0,&RAP)); 389*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMatMult(R,A,P,MAT_REUSE_MATRIX,2.0,&RAP)); 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNormDifference(C,RAP,&norm)); 3912c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PtAP != RAP %g",(double)norm); 392*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&R)); 393*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RAP)); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown /* Create vector x that is compatible with P */ 397*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(P,&m,&n)); 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,n,PETSC_DECIDE)); 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 401c4762a1bSJed Brown 402*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&v3)); 403*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(v3,n,PETSC_DECIDE)); 404*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(v3)); 405*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(v3,&v4)); 406c4762a1bSJed Brown 407c4762a1bSJed Brown norm = 0.0; 408c4762a1bSJed Brown for (i=0; i<10; i++) { 409*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 410*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(P,x,v1)); 411*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,v1,v2)); /* v2 = A*P*x */ 412c4762a1bSJed Brown 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(P,v2,v3)); /* v3 = Pt*A*P*x */ 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,v4)); /* v3 = C*x */ 415*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v4,NORM_2,&norm_abs)); 416*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v4,none,v3)); 417*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(v4,NORM_2,&norm_tmp)); 418c4762a1bSJed Brown 419c4762a1bSJed Brown norm_tmp /= norm_abs; 420c4762a1bSJed Brown if (norm_tmp > norm) norm = norm_tmp; 421c4762a1bSJed Brown } 4222c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm >= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error: MatPtAP(), |v1 - v2|: %g",(double)norm); 423c4762a1bSJed Brown 424*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 425*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 426*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v3)); 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v4)); 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* Destroy objects */ 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v1)); 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v2)); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 436*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idxn)); 437c4762a1bSJed Brown 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A_save)); 439*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 440c4762a1bSJed Brown 441c4762a1bSJed Brown PetscPreLoadEnd(); 442c4762a1bSJed Brown PetscFinalize(); 443c4762a1bSJed Brown return ierr; 444c4762a1bSJed Brown } 445c4762a1bSJed Brown 446c4762a1bSJed Brown /*TEST 447c4762a1bSJed Brown 448c4762a1bSJed Brown test: 449c4762a1bSJed Brown suffix: 2_mattransposematmult_matmatmult 450c4762a1bSJed Brown nsize: 3 451dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 452c20d7725SJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via at*b> ex94_2.tmp 2>&1 453c4762a1bSJed Brown 454c4762a1bSJed Brown test: 455c4762a1bSJed Brown suffix: 2_mattransposematmult_scalable 456c4762a1bSJed Brown nsize: 3 457dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 458c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -mattransposematmult_via scalable> ex94_2.tmp 2>&1 459c4762a1bSJed Brown output_file: output/ex94_1.out 460c4762a1bSJed Brown 461c4762a1bSJed Brown test: 462c4762a1bSJed Brown suffix: axpy_mpiaij 463dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 464c4762a1bSJed Brown nsize: 8 465c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY 466c4762a1bSJed Brown output_file: output/ex94_1.out 467c4762a1bSJed Brown 468c4762a1bSJed Brown test: 469c4762a1bSJed Brown suffix: axpy_mpibaij 470dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 471c4762a1bSJed Brown nsize: 8 472c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type baij 473c4762a1bSJed Brown output_file: output/ex94_1.out 474c4762a1bSJed Brown 475c4762a1bSJed Brown test: 476c4762a1bSJed Brown suffix: axpy_mpisbaij 477dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 478c4762a1bSJed Brown nsize: 8 479c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/poisson_2d5p -f1 ${DATAFILESPATH}/matrices/poisson_2d13p -test_MatAXPY -mat_type sbaij 480c4762a1bSJed Brown output_file: output/ex94_1.out 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: matmatmult 484dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 485c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info 486c4762a1bSJed Brown output_file: output/ex94_1.out 487c4762a1bSJed Brown 488c4762a1bSJed Brown test: 489c4762a1bSJed Brown suffix: matmatmult_2 490dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 491c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -mat_type mpiaij -viewer_binary_skip_info 492c4762a1bSJed Brown output_file: output/ex94_1.out 493c4762a1bSJed Brown 494c4762a1bSJed Brown test: 495c4762a1bSJed Brown suffix: matmatmult_scalable 496c4762a1bSJed Brown nsize: 4 497dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 498c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -matmatmult_via scalable 499c4762a1bSJed Brown output_file: output/ex94_1.out 500c4762a1bSJed Brown 501c4762a1bSJed Brown test: 502c4762a1bSJed Brown suffix: ptap 503c4762a1bSJed Brown nsize: 3 504dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 505c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium -matptap_via scalable 506c4762a1bSJed Brown output_file: output/ex94_1.out 507c4762a1bSJed Brown 508c4762a1bSJed Brown test: 509c4762a1bSJed Brown suffix: rap 510c4762a1bSJed Brown nsize: 3 511dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 512c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/medium -f1 ${DATAFILESPATH}/matrices/medium 513c4762a1bSJed Brown output_file: output/ex94_1.out 514c4762a1bSJed Brown 515c4762a1bSJed Brown test: 516c4762a1bSJed Brown suffix: scalable0 517dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 518c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info 519c4762a1bSJed Brown output_file: output/ex94_1.out 520c4762a1bSJed Brown 521c4762a1bSJed Brown test: 522c4762a1bSJed Brown suffix: scalable1 523dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 524c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/arco1 -f1 ${DATAFILESPATH}/matrices/arco1 -viewer_binary_skip_info -matptap_via scalable 525c4762a1bSJed Brown output_file: output/ex94_1.out 526c4762a1bSJed Brown 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: view 529c4762a1bSJed Brown nsize: 2 530dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 531c4762a1bSJed Brown args: -f0 ${DATAFILESPATH}/matrices/tiny -f1 ${DATAFILESPATH}/matrices/tiny -viewer_binary_skip_info -matops_view 532c4762a1bSJed Brown output_file: output/ex94_2.out 533c4762a1bSJed Brown 534c4762a1bSJed Brown TEST*/ 535