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