1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatTransposeMatMult() on MatLoad() matrix \n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **args) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat A,C,Bdense,Cdense; 9*c4762a1bSJed Brown PetscErrorCode ierr; 10*c4762a1bSJed Brown PetscViewer fd; /* viewer */ 11*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 12*c4762a1bSJed Brown PetscBool flg,viewmats=PETSC_FALSE; 13*c4762a1bSJed Brown PetscMPIInt rank,size; 14*c4762a1bSJed Brown PetscReal fill=1.0; 15*c4762a1bSJed Brown PetscInt m,n,i,j,BN=10,rstart,rend,*rows,*cols; 16*c4762a1bSJed Brown PetscScalar *Barray,*Carray,rval,*array; 17*c4762a1bSJed Brown Vec x,y; 18*c4762a1bSJed Brown PetscRandom rand; 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 21*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 24*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 25*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown /* Load matrix A */ 28*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown /* Print (for testing only) */ 34*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats);CHKERRQ(ierr); 35*c4762a1bSJed Brown if (viewmats) { 36*c4762a1bSJed Brown if (!rank) printf("A_aij:\n"); 37*c4762a1bSJed Brown ierr = MatView(A,0);CHKERRQ(ierr); 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_aij() */ 41*c4762a1bSJed Brown ierr = MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 42*c4762a1bSJed Brown if (viewmats) { 43*c4762a1bSJed Brown if (!rank) printf("\nC = A_aij^T * A_aij:\n"); 44*c4762a1bSJed Brown ierr = MatView(C,0);CHKERRQ(ierr); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown /* create a dense matrix Bdense */ 50*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Bdense);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = MatSetType(Bdense,MATDENSE);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatSetFromOptions(Bdense);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatSetUp(Bdense);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatGetOwnershipRange(Bdense,&rstart,&rend);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown ierr = PetscMalloc3(m,&rows,BN,&cols,m*BN,&array);CHKERRQ(ierr); 58*c4762a1bSJed Brown for (i=0; i<m; i++) rows[i] = rstart + i; 59*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 61*c4762a1bSJed Brown for (j=0; j<BN; j++) { 62*c4762a1bSJed Brown cols[j] = j; 63*c4762a1bSJed Brown for (i=0; i<m; i++) { 64*c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 65*c4762a1bSJed Brown array[m*j+i] = rval; 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown ierr = MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscFree3(rows,cols,array);CHKERRQ(ierr); 73*c4762a1bSJed Brown if (viewmats) { 74*c4762a1bSJed Brown if (!rank) printf("\nBdense:\n"); 75*c4762a1bSJed Brown ierr = MatView(Bdense,0);CHKERRQ(ierr); 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown /* Test MatTransposeMatMult_aij_dense() */ 79*c4762a1bSJed Brown ierr = MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 81*c4762a1bSJed Brown if (viewmats) { 82*c4762a1bSJed Brown if (!rank) printf("\nC=A^T*Bdense:\n"); 83*c4762a1bSJed Brown ierr = MatView(C,0);CHKERRQ(ierr); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* Check accuracy */ 87*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Cdense);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatSetType(Cdense,MATDENSE);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatSetFromOptions(Cdense);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatSetUp(Cdense);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 96*c4762a1bSJed Brown if (size == 1) { 97*c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y);CHKERRQ(ierr); 99*c4762a1bSJed Brown } else { 100*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y);CHKERRQ(ierr); 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown /* Cdense[:,j] = A^T * Bdense[:,j] */ 105*c4762a1bSJed Brown ierr = MatDenseGetArray(Bdense,&Barray);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = MatDenseGetArray(Cdense,&Carray);CHKERRQ(ierr); 107*c4762a1bSJed Brown for (j=0; j<BN; j++) { 108*c4762a1bSJed Brown ierr = VecPlaceArray(x,Barray);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecPlaceArray(y,Carray);CHKERRQ(ierr); 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown ierr = MatMultTranspose(A,x,y);CHKERRQ(ierr); 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown ierr = VecResetArray(x);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = VecResetArray(y);CHKERRQ(ierr); 115*c4762a1bSJed Brown Barray += m; 116*c4762a1bSJed Brown Carray += n; 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown ierr = MatDenseRestoreArray(Bdense,&Barray);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = MatDenseRestoreArray(Cdense,&Carray);CHKERRQ(ierr); 120*c4762a1bSJed Brown if (viewmats) { 121*c4762a1bSJed Brown if (!rank) printf("\nCdense:\n"); 122*c4762a1bSJed Brown ierr = MatView(Cdense,0);CHKERRQ(ierr); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = MatEqual(C,Cdense,&flg);CHKERRQ(ierr); 126*c4762a1bSJed Brown if (!flg) { 127*c4762a1bSJed Brown if (!rank) printf(" C != Cdense\n"); 128*c4762a1bSJed Brown } 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown /* Free data structures */ 131*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = MatDestroy(&Bdense);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = MatDestroy(&Cdense);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = PetscFinalize(); 138*c4762a1bSJed Brown return ierr; 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /*TEST 143*c4762a1bSJed Brown 144*c4762a1bSJed Brown test: 145*c4762a1bSJed Brown requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES) 146*c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 147*c4762a1bSJed Brown output_file: output/ex163.out 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown test: 150*c4762a1bSJed Brown suffix: 2 151*c4762a1bSJed Brown nsize: 3 152*c4762a1bSJed Brown requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES) 153*c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small 154*c4762a1bSJed Brown output_file: output/ex163.out 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown TEST*/ 157