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; 21ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 24589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 25c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* Load matrix A */ 28c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Print (for testing only) */ 34c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats);CHKERRQ(ierr); 35c4762a1bSJed Brown if (viewmats) { 36c4762a1bSJed Brown if (!rank) printf("A_aij:\n"); 37c4762a1bSJed Brown ierr = MatView(A,0);CHKERRQ(ierr); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_aij() */ 41c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 42c4762a1bSJed Brown if (viewmats) { 43c4762a1bSJed Brown if (!rank) printf("\nC = A_aij^T * A_aij:\n"); 44c4762a1bSJed Brown ierr = MatView(C,0);CHKERRQ(ierr); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* create a dense matrix Bdense */ 50c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Bdense);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatSetType(Bdense,MATDENSE);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = MatSetFromOptions(Bdense);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatSetUp(Bdense);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = MatGetOwnershipRange(Bdense,&rstart,&rend);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscMalloc3(m,&rows,BN,&cols,m*BN,&array);CHKERRQ(ierr); 58c4762a1bSJed Brown for (i=0; i<m; i++) rows[i] = rstart + i; 59c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 61c4762a1bSJed Brown for (j=0; j<BN; j++) { 62c4762a1bSJed Brown cols[j] = j; 63c4762a1bSJed Brown for (i=0; i<m; i++) { 64c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 65c4762a1bSJed Brown array[m*j+i] = rval; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown } 68c4762a1bSJed Brown ierr = MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = PetscFree3(rows,cols,array);CHKERRQ(ierr); 73c4762a1bSJed Brown if (viewmats) { 74c4762a1bSJed Brown if (!rank) printf("\nBdense:\n"); 75c4762a1bSJed Brown ierr = MatView(Bdense,0);CHKERRQ(ierr); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_dense() */ 79c4762a1bSJed Brown ierr = MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 81c4762a1bSJed Brown if (viewmats) { 82c4762a1bSJed Brown if (!rank) printf("\nC=A^T*Bdense:\n"); 83c4762a1bSJed Brown ierr = MatView(C,0);CHKERRQ(ierr); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Check accuracy */ 87c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Cdense);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = MatSetType(Cdense,MATDENSE);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatSetFromOptions(Cdense);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatSetUp(Cdense);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94c4762a1bSJed Brown 95ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 96c4762a1bSJed Brown if (size == 1) { 97c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y);CHKERRQ(ierr); 99c4762a1bSJed Brown } else { 100c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y);CHKERRQ(ierr); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Cdense[:,j] = A^T * Bdense[:,j] */ 105c4762a1bSJed Brown ierr = MatDenseGetArray(Bdense,&Barray);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatDenseGetArray(Cdense,&Carray);CHKERRQ(ierr); 107c4762a1bSJed Brown for (j=0; j<BN; j++) { 108c4762a1bSJed Brown ierr = VecPlaceArray(x,Barray);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = VecPlaceArray(y,Carray);CHKERRQ(ierr); 110c4762a1bSJed Brown 111c4762a1bSJed Brown ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown ierr = VecResetArray(x);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = VecResetArray(y);CHKERRQ(ierr); 115c4762a1bSJed Brown Barray += m; 116c4762a1bSJed Brown Carray += n; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown ierr = MatDenseRestoreArray(Bdense,&Barray);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDenseRestoreArray(Cdense,&Carray);CHKERRQ(ierr); 120c4762a1bSJed Brown if (viewmats) { 121c4762a1bSJed Brown if (!rank) printf("\nCdense:\n"); 122c4762a1bSJed Brown ierr = MatView(Cdense,0);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown ierr = MatEqual(C,Cdense,&flg);CHKERRQ(ierr); 126c4762a1bSJed Brown if (!flg) { 127c4762a1bSJed Brown if (!rank) printf(" C != Cdense\n"); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Free data structures */ 131c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = MatDestroy(&Bdense);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = MatDestroy(&Cdense);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = PetscFinalize(); 138c4762a1bSJed Brown return ierr; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /*TEST 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144*dfd57a17SPierre 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 151*dfd57a17SPierre 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