1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **argv) 6c4762a1bSJed Brown { 7c20d7725SJed Brown Mat A,B,C,D,AT; 8c4762a1bSJed Brown PetscInt i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an; 9c4762a1bSJed Brown PetscScalar v; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscRandom r; 12*3fbca975SStefano Zampini PetscBool equal=PETSC_FALSE,flg; 13c20d7725SJed Brown PetscReal fill = 1.0,norm; 14c4762a1bSJed Brown PetscMPIInt size; 15c4762a1bSJed Brown 16c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 17c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 18c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 19c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL);CHKERRQ(ierr); 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 22c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* Create a aij matrix A */ 25c4762a1bSJed Brown M = N = m*n; 26c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,5,NULL,5,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,5,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr); 34c4762a1bSJed Brown am = Iend - Istart; 35c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 36c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 37c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 38c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 39c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 40c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 41c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Create a dense matrix B */ 47c4762a1bSJed Brown ierr = MatGetLocalSize(A,&am,&an);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = MatMPIDenseSetPreallocation(B,NULL);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = MatSetRandom(B,r);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 58c4762a1bSJed Brown 59c20d7725SJed Brown /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */ 60c20d7725SJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 61c20d7725SJed Brown ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr); 62c20d7725SJed Brown ierr = MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N);CHKERRQ(ierr); 63c20d7725SJed Brown ierr = MatSetFromOptions(C); CHKERRQ(ierr); 64c20d7725SJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 65c20d7725SJed Brown ierr = MatZeroEntries(C);CHKERRQ(ierr); 66c20d7725SJed Brown ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 67c20d7725SJed Brown ierr = MatNorm(C,NORM_INFINITY,&norm);CHKERRQ(ierr); 68c20d7725SJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 69c20d7725SJed Brown 70c4762a1bSJed Brown /* Test C = A*B (aij*dense) */ 71c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); 73c4762a1bSJed Brown 74c20d7725SJed Brown /* Test developer API */ 75c20d7725SJed Brown ierr = MatProductCreate(A,B,NULL,&D);CHKERRQ(ierr); 76c20d7725SJed Brown ierr = MatProductSetType(D,MATPRODUCT_AB);CHKERRQ(ierr); 77c20d7725SJed Brown ierr = MatProductSetAlgorithm(D,"default");CHKERRQ(ierr); 78c20d7725SJed Brown ierr = MatProductSetFill(D,fill);CHKERRQ(ierr); 79c20d7725SJed Brown ierr = MatProductSetFromOptions(D);CHKERRQ(ierr); 80c20d7725SJed Brown ierr = MatProductSymbolic(D);CHKERRQ(ierr); 81c4762a1bSJed Brown for (i=0; i<2; i++) { 82c20d7725SJed Brown ierr = MatProductNumeric(D);CHKERRQ(ierr); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown ierr = MatEqual(C,D,&equal);CHKERRQ(ierr); 85c4762a1bSJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D"); 86c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 87c4762a1bSJed Brown 88c20d7725SJed Brown /* Test D = AT*B (transpose(aij)*dense) */ 89c20d7725SJed Brown ierr = MatCreateTranspose(A,&AT);CHKERRQ(ierr); 90c20d7725SJed Brown ierr = MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 91c20d7725SJed Brown ierr = MatMatMultEqual(AT,B,D,10,&equal);CHKERRQ(ierr); 92c20d7725SJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)"); 93c20d7725SJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 94c20d7725SJed Brown ierr = MatDestroy(&AT);CHKERRQ(ierr); 95c20d7725SJed Brown 96c4762a1bSJed Brown /* Test D = C*A (dense*aij) */ 97c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 99c20d7725SJed Brown ierr = MatMatMultEqual(C,A,D,10,&equal);CHKERRQ(ierr); 100c20d7725SJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)"); 101c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test D = A*C (aij*dense) */ 104c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 106c20d7725SJed Brown ierr = MatMatMultEqual(A,C,D,10,&equal);CHKERRQ(ierr); 107c20d7725SJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)"); 108c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Test D = B*C (dense*dense) */ 111c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 112c4762a1bSJed Brown if (size == 1) { 113c4762a1bSJed Brown ierr = MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 115c20d7725SJed Brown ierr = MatMatMultEqual(B,C,D,10,&equal);CHKERRQ(ierr); 116c20d7725SJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)"); 117c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c20d7725SJed Brown /* Test D = B*C^T (dense*dense) */ 121c20d7725SJed Brown ierr = MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 122c20d7725SJed Brown ierr = MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D);CHKERRQ(ierr); 123c20d7725SJed Brown ierr = MatMatTransposeMultEqual(B,C,D,10,&equal);CHKERRQ(ierr); 124c20d7725SJed Brown if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)"); 125c20d7725SJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 126c20d7725SJed Brown 127*3fbca975SStefano Zampini /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */ 128*3fbca975SStefano Zampini flg = PETSC_FALSE; 129*3fbca975SStefano Zampini ierr = PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg);CHKERRQ(ierr); 130*3fbca975SStefano Zampini if (flg) { 131*3fbca975SStefano Zampini /* user driver */ 132*3fbca975SStefano Zampini ierr = MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B);CHKERRQ(ierr); 133*3fbca975SStefano Zampini } else { 134*3fbca975SStefano Zampini /* clear internal data structures related with previous products to avoid circular references */ 135*3fbca975SStefano Zampini ierr = MatProductClear(A);CHKERRQ(ierr); 136*3fbca975SStefano Zampini ierr = MatProductClear(B);CHKERRQ(ierr); 137*3fbca975SStefano Zampini ierr = MatProductClear(C);CHKERRQ(ierr); 138*3fbca975SStefano Zampini ierr = MatProductCreateWithMat(A,C,NULL,B);CHKERRQ(ierr); 139*3fbca975SStefano Zampini ierr = MatProductSetType(B,MATPRODUCT_AB);CHKERRQ(ierr); 140*3fbca975SStefano Zampini ierr = MatProductSetFromOptions(B);CHKERRQ(ierr); 141*3fbca975SStefano Zampini ierr = MatProductNumeric(B);CHKERRQ(ierr); 142*3fbca975SStefano Zampini } 143*3fbca975SStefano Zampini 144c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = PetscFinalize(); 148c4762a1bSJed Brown return ierr; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown 152c4762a1bSJed Brown /*TEST 153c4762a1bSJed Brown 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown args: -M 10 -N 10 156c4762a1bSJed Brown output_file: output/ex109.out 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: 2 160c4762a1bSJed Brown nsize: 3 161c4762a1bSJed Brown output_file: output/ex109.out 162c4762a1bSJed Brown 163c20d7725SJed Brown test: 164c20d7725SJed Brown suffix: 3 165c20d7725SJed Brown nsize: 2 166c20d7725SJed Brown args: -matmattransmult_mpidense_mpidense_via cyclic 167c20d7725SJed Brown output_file: output/ex109.out 168c20d7725SJed Brown 169*3fbca975SStefano Zampini test: 170*3fbca975SStefano Zampini suffix: 4 171*3fbca975SStefano Zampini args: -test_userAPI 172*3fbca975SStefano Zampini output_file: output/ex109.out 173*3fbca975SStefano Zampini 174*3fbca975SStefano Zampini test: 175*3fbca975SStefano Zampini suffix: 5 176*3fbca975SStefano Zampini nsize: 3 177*3fbca975SStefano Zampini args: -test_userAPI 178*3fbca975SStefano Zampini output_file: output/ex109.out 179*3fbca975SStefano Zampini 180c4762a1bSJed Brown TEST*/ 181