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 PetscErrorCode ierr; 10c4762a1bSJed Brown PetscViewer fd; /* viewer */ 11c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 12c4762a1bSJed Brown PetscBool flg,viewmats=PETSC_FALSE; 13c4762a1bSJed Brown PetscMPIInt rank,size; 14c4762a1bSJed Brown PetscReal fill=1.0; 15c4762a1bSJed Brown PetscInt m,n,i,j,BN=10,rstart,rend,*rows,*cols; 16c4762a1bSJed Brown PetscScalar *Barray,*Carray,rval,*array; 17c4762a1bSJed Brown Vec x,y; 18c4762a1bSJed Brown PetscRandom rand; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 25*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* Load matrix A */ 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Print (for testing only) */ 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats)); 35c4762a1bSJed Brown if (viewmats) { 36dd400576SPatrick Sanan if (rank == 0) printf("A_aij:\n"); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,0)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_aij() */ 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C)); 42c4762a1bSJed Brown if (viewmats) { 43dd400576SPatrick Sanan if (rank == 0) printf("\nC = A_aij^T * A_aij:\n"); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,0)); 45c4762a1bSJed Brown } 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* create a dense matrix Bdense */ 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Bdense)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Bdense,MATDENSE)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Bdense)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Bdense)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(Bdense,&rstart,&rend)); 56c4762a1bSJed Brown 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(m,&rows,BN,&cols,m*BN,&array)); 58c4762a1bSJed Brown for (i=0; i<m; i++) rows[i] = rstart + i; 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 61c4762a1bSJed Brown for (j=0; j<BN; j++) { 62c4762a1bSJed Brown cols[j] = j; 63c4762a1bSJed Brown for (i=0; i<m; i++) { 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 65c4762a1bSJed Brown array[m*j+i] = rval; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown } 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(rows,cols,array)); 73c4762a1bSJed Brown if (viewmats) { 74dd400576SPatrick Sanan if (rank == 0) printf("\nBdense:\n"); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Bdense,0)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_dense() */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C)); 81c4762a1bSJed Brown if (viewmats) { 82dd400576SPatrick Sanan if (rank == 0) printf("\nC=A^T*Bdense:\n"); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,0)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Check accuracy */ 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Cdense)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Cdense,MATDENSE)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Cdense)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Cdense)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 96c4762a1bSJed Brown if (size == 1) { 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y)); 99c4762a1bSJed Brown } else { 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Cdense[:,j] = A^T * Bdense[:,j] */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(Bdense,&Barray)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(Cdense,&Carray)); 107c4762a1bSJed Brown for (j=0; j<BN; j++) { 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(x,Barray)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecPlaceArray(y,Carray)); 110c4762a1bSJed Brown 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,x,y)); 112c4762a1bSJed Brown 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(x)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecResetArray(y)); 115c4762a1bSJed Brown Barray += m; 116c4762a1bSJed Brown Carray += n; 117c4762a1bSJed Brown } 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(Bdense,&Barray)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(Cdense,&Carray)); 120c4762a1bSJed Brown if (viewmats) { 121dd400576SPatrick Sanan if (rank == 0) printf("\nCdense:\n"); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Cdense,0)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(C,Cdense,&flg)); 126c4762a1bSJed Brown if (!flg) { 127dd400576SPatrick Sanan if (rank == 0) printf(" C != Cdense\n"); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Free data structures */ 1315f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Bdense)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 137c4762a1bSJed Brown ierr = PetscFinalize(); 138c4762a1bSJed Brown return ierr; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /*TEST 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 145c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 146c4762a1bSJed Brown output_file: output/ex163.out 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: 2 150c4762a1bSJed Brown nsize: 3 151dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 152c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 153c4762a1bSJed Brown output_file: output/ex163.out 154c4762a1bSJed Brown 155c4762a1bSJed Brown TEST*/ 156