1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatTransposeMatMult() on MatLoad() matrix \n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,C,Bdense,Cdense; 9c4762a1bSJed Brown PetscViewer fd; /* viewer */ 10c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 11c4762a1bSJed Brown PetscBool flg,viewmats=PETSC_FALSE; 12c4762a1bSJed Brown PetscMPIInt rank,size; 13c4762a1bSJed Brown PetscReal fill=1.0; 14c4762a1bSJed Brown PetscInt m,n,i,j,BN=10,rstart,rend,*rows,*cols; 15c4762a1bSJed Brown PetscScalar *Barray,*Carray,rval,*array; 16c4762a1bSJed Brown Vec x,y; 17c4762a1bSJed Brown PetscRandom rand; 18c4762a1bSJed Brown 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 2428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Load matrix A */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Print (for testing only) */ 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats)); 34c4762a1bSJed Brown if (viewmats) { 35dd400576SPatrick Sanan if (rank == 0) printf("A_aij:\n"); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,0)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_aij() */ 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C)); 41c4762a1bSJed Brown if (viewmats) { 42dd400576SPatrick Sanan if (rank == 0) printf("\nC = A_aij^T * A_aij:\n"); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,0)); 44c4762a1bSJed Brown } 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* create a dense matrix Bdense */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Bdense)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Bdense,MATDENSE)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Bdense)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Bdense)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(Bdense,&rstart,&rend)); 55c4762a1bSJed Brown 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(m,&rows,BN,&cols,m*BN,&array)); 57c4762a1bSJed Brown for (i=0; i<m; i++) rows[i] = rstart + i; 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 60c4762a1bSJed Brown for (j=0; j<BN; j++) { 61c4762a1bSJed Brown cols[j] = j; 62c4762a1bSJed Brown for (i=0; i<m; i++) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 64c4762a1bSJed Brown array[m*j+i] = rval; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown } 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(rows,cols,array)); 72c4762a1bSJed Brown if (viewmats) { 73dd400576SPatrick Sanan if (rank == 0) printf("\nBdense:\n"); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Bdense,0)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_dense() */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C)); 80c4762a1bSJed Brown if (viewmats) { 81dd400576SPatrick Sanan if (rank == 0) printf("\nC=A^T*Bdense:\n"); 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,0)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Check accuracy */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Cdense)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Cdense,MATDENSE)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Cdense)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Cdense)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY)); 93c4762a1bSJed Brown 945f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 95c4762a1bSJed Brown if (size == 1) { 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y)); 98c4762a1bSJed Brown } else { 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Cdense[:,j] = A^T * Bdense[:,j] */ 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(Bdense,&Barray)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(Cdense,&Carray)); 106c4762a1bSJed Brown for (j=0; j<BN; j++) { 1075f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(x,Barray)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(y,Carray)); 109c4762a1bSJed Brown 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,y)); 111c4762a1bSJed Brown 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(x)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(y)); 114c4762a1bSJed Brown Barray += m; 115c4762a1bSJed Brown Carray += n; 116c4762a1bSJed Brown } 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(Bdense,&Barray)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(Cdense,&Carray)); 119c4762a1bSJed Brown if (viewmats) { 120dd400576SPatrick Sanan if (rank == 0) printf("\nCdense:\n"); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Cdense,0)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(C,Cdense,&flg)); 125c4762a1bSJed Brown if (!flg) { 126dd400576SPatrick Sanan if (rank == 0) printf(" C != Cdense\n"); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Free data structures */ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Bdense)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 136*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 137*b122ec5aSJacob Faibussowitsch return 0; 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown /*TEST 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 144c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 145c4762a1bSJed Brown output_file: output/ex163.out 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: 2 149c4762a1bSJed Brown nsize: 3 150dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 151c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 152c4762a1bSJed Brown output_file: output/ex163.out 153c4762a1bSJed Brown 154c4762a1bSJed Brown TEST*/ 155