1076ba34aSJunchao Zhang static char help[] = "Test Mat products \n\n"; 2076ba34aSJunchao Zhang 3076ba34aSJunchao Zhang #include <petscmat.h> 4076ba34aSJunchao Zhang int main(int argc,char **args) 5076ba34aSJunchao Zhang { 6076ba34aSJunchao Zhang Mat A=NULL,B=NULL,C=NULL,D=NULL,E=NULL; 7076ba34aSJunchao Zhang PetscErrorCode ierr; 8076ba34aSJunchao Zhang PetscInt k; 9076ba34aSJunchao Zhang const PetscInt M = 18,N = 18; 10076ba34aSJunchao Zhang PetscMPIInt rank; 11076ba34aSJunchao Zhang 12076ba34aSJunchao Zhang /* A, B are 18 x 18 nonsymmetric matrices and have the same sparsity pattern but different values. 13076ba34aSJunchao Zhang Big enough to have complex communication patterns but still small enough for debugging. 14076ba34aSJunchao Zhang */ 15076ba34aSJunchao Zhang 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}; 16076ba34aSJunchao Zhang 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}; 17076ba34aSJunchao Zhang 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}; 18076ba34aSJunchao Zhang 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}; 19076ba34aSJunchao Zhang 20076ba34aSJunchao Zhang PetscInt Annz = sizeof(Ai)/sizeof(PetscInt); 21076ba34aSJunchao Zhang PetscInt Bnnz = sizeof(Bi)/sizeof(PetscInt); 22076ba34aSJunchao Zhang 23076ba34aSJunchao Zhang ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 25076ba34aSJunchao Zhang 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 27076ba34aSJunchao Zhang ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 29076ba34aSJunchao Zhang ierr = MatSeqAIJSetPreallocation(A,2,NULL); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,2,NULL,2,NULL)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 32076ba34aSJunchao Zhang 33076ba34aSJunchao Zhang if (rank == 0) { 34*5f80ce2aSJacob Faibussowitsch for (k=0; k<Annz; k++) CHKERRQ(MatSetValue(A,Ai[k],Aj[k],Ai[k]+Aj[k]+1.0,INSERT_VALUES)); 35076ba34aSJunchao Zhang } 36076ba34aSJunchao Zhang 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 39076ba34aSJunchao Zhang 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 41076ba34aSJunchao Zhang ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 43076ba34aSJunchao Zhang ierr = MatSeqAIJSetPreallocation(B,2,NULL); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(B,2,NULL,2,NULL)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 46076ba34aSJunchao Zhang 47076ba34aSJunchao Zhang if (rank == 0) { 48*5f80ce2aSJacob Faibussowitsch for (k=0; k<Bnnz; k++) CHKERRQ(MatSetValue(B,Bi[k],Bj[k],Bi[k]+Bj[k]+2.0,INSERT_VALUES)); 49076ba34aSJunchao Zhang } 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 52076ba34aSJunchao Zhang 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 55076ba34aSJunchao Zhang 56076ba34aSJunchao Zhang /* B, A have the same nonzero pattern, so it is legitimate to do so */ 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 59076ba34aSJunchao Zhang 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(D, PETSC_VIEWER_STDOUT_WORLD)); 62076ba34aSJunchao Zhang 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&E)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(E,PETSC_VIEWER_STDOUT_WORLD)); 65076ba34aSJunchao Zhang 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&E)); 71076ba34aSJunchao Zhang 72076ba34aSJunchao Zhang ierr = PetscFinalize(); 73076ba34aSJunchao Zhang return ierr; 74076ba34aSJunchao Zhang } 75076ba34aSJunchao Zhang 76076ba34aSJunchao Zhang /*TEST 77076ba34aSJunchao Zhang testset: 78076ba34aSJunchao Zhang filter: grep -ve type -ve "Mat Object" 79076ba34aSJunchao Zhang output_file: output/ex250_1.out 80076ba34aSJunchao Zhang 81076ba34aSJunchao Zhang test: 82076ba34aSJunchao Zhang suffix: 1 83076ba34aSJunchao Zhang nsize: {{1 3}} 84076ba34aSJunchao Zhang args: -mat_type aij 85076ba34aSJunchao Zhang 86076ba34aSJunchao Zhang test: 87076ba34aSJunchao Zhang suffix: 2 88076ba34aSJunchao Zhang nsize: {{3 4}} 89076ba34aSJunchao Zhang args: -mat_type aij -matmatmult_via backend -matptap_via backend -mattransposematmult_via backend 90076ba34aSJunchao Zhang 91076ba34aSJunchao Zhang test: 92076ba34aSJunchao Zhang suffix: cuda 93076ba34aSJunchao Zhang requires: cuda 94076ba34aSJunchao Zhang nsize: {{1 3 4}} 95076ba34aSJunchao Zhang args: -mat_type aijcusparse 96076ba34aSJunchao Zhang 97076ba34aSJunchao Zhang test: 98076ba34aSJunchao Zhang suffix: kok 993078479eSJunchao Zhang requires: !sycl kokkos_kernels 100076ba34aSJunchao Zhang nsize: {{1 3 4}} 101076ba34aSJunchao Zhang args: -mat_type aijkokkos 102076ba34aSJunchao Zhang 103076ba34aSJunchao Zhang TEST*/ 104