1 static char help[] = "Test Mat products \n\n"; 2 3 #include <petscmat.h> 4 int main(int argc,char **args) 5 { 6 Mat A=NULL,B=NULL,C=NULL,D=NULL,E=NULL; 7 PetscErrorCode ierr; 8 PetscInt k; 9 const PetscInt M = 18,N = 18; 10 PetscMPIInt rank; 11 12 /* A, B are 18 x 18 nonsymmetric matrices and have the same sparsity pattern but different values. 13 Big enough to have complex communication patterns but still small enough for debugging. 14 */ 15 PetscInt Ai[] = {0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17}; 16 PetscInt Aj[] = {0, 1, 2, 7, 3, 8, 4, 9, 5, 8, 2, 6, 11, 0, 7, 1, 6, 2, 4, 10, 16, 11, 15, 12, 17, 12, 13, 14, 15, 17, 11, 13, 3, 16, 9, 15, 11, 13}; 17 PetscInt Bi[] = {0, 0, 1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11, 11, 12, 12, 13, 13, 13, 14, 14, 15, 15, 16, 16, 17, 17}; 18 PetscInt Bj[] = {0, 1, 2, 7, 3, 8, 4, 9, 5, 8, 2, 6, 11, 0, 7, 1, 6, 2, 4, 10, 16, 11, 15, 12, 17, 12, 13, 14, 15, 17, 11, 13, 3, 16, 9, 15, 11, 13}; 19 20 PetscInt Annz = sizeof(Ai)/sizeof(PetscInt); 21 PetscInt Bnnz = sizeof(Bi)/sizeof(PetscInt); 22 23 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 25 26 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 27 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N); 28 CHKERRQ(MatSetFromOptions(A)); 29 ierr = MatSeqAIJSetPreallocation(A,2,NULL); 30 CHKERRQ(MatMPIAIJSetPreallocation(A,2,NULL,2,NULL)); 31 CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 32 33 if (rank == 0) { 34 for (k=0; k<Annz; k++) CHKERRQ(MatSetValue(A,Ai[k],Aj[k],Ai[k]+Aj[k]+1.0,INSERT_VALUES)); 35 } 36 37 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 38 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 39 40 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 41 ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N); 42 CHKERRQ(MatSetFromOptions(B)); 43 ierr = MatSeqAIJSetPreallocation(B,2,NULL); 44 CHKERRQ(MatMPIAIJSetPreallocation(B,2,NULL,2,NULL)); 45 CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 46 47 if (rank == 0) { 48 for (k=0; k<Bnnz; k++) CHKERRQ(MatSetValue(B,Bi[k],Bj[k],Bi[k]+Bj[k]+2.0,INSERT_VALUES)); 49 } 50 CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 51 CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 52 53 CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 54 CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 55 56 /* B, A have the same nonzero pattern, so it is legitimate to do so */ 57 CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 58 CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 59 60 CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D)); 61 CHKERRQ(MatView(D, PETSC_VIEWER_STDOUT_WORLD)); 62 63 CHKERRQ(MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&E)); 64 CHKERRQ(MatView(E,PETSC_VIEWER_STDOUT_WORLD)); 65 66 CHKERRQ(MatDestroy(&A)); 67 CHKERRQ(MatDestroy(&B)); 68 CHKERRQ(MatDestroy(&C)); 69 CHKERRQ(MatDestroy(&D)); 70 CHKERRQ(MatDestroy(&E)); 71 72 ierr = PetscFinalize(); 73 return ierr; 74 } 75 76 /*TEST 77 testset: 78 filter: grep -ve type -ve "Mat Object" 79 output_file: output/ex250_1.out 80 81 test: 82 suffix: 1 83 nsize: {{1 3}} 84 args: -mat_type aij 85 86 test: 87 suffix: 2 88 nsize: {{3 4}} 89 args: -mat_type aij -matmatmult_via backend -matptap_via backend -mattransposematmult_via backend 90 91 test: 92 suffix: cuda 93 requires: cuda 94 nsize: {{1 3 4}} 95 args: -mat_type aijcusparse 96 97 test: 98 suffix: kok 99 requires: !sycl kokkos_kernels 100 nsize: {{1 3 4}} 101 args: -mat_type aijkokkos 102 103 TEST*/ 104