121c60b78SJunchao Zhang static char help[] = "Test MatSetValuesCOO for MPIAIJ and its subclasses \n\n"; 242550becSJunchao Zhang 342550becSJunchao Zhang #include <petscmat.h> 442550becSJunchao Zhang int main(int argc,char **args) 542550becSJunchao Zhang { 642550becSJunchao Zhang Mat A,B; 742550becSJunchao Zhang PetscInt k; 842550becSJunchao Zhang const PetscInt M = 18,N = 18; 942550becSJunchao Zhang PetscMPIInt rank,size; 1042550becSJunchao Zhang PetscBool equal; 1142550becSJunchao Zhang PetscScalar *vals; 1221c60b78SJunchao Zhang PetscBool flg = PETSC_FALSE; 1342550becSJunchao Zhang 1442550becSJunchao Zhang /* Construct 18 x 18 matrices, which are big enough to have complex communication patterns but still small enough for debugging */ 1542550becSJunchao Zhang PetscInt i0[] = {7, 7, 8, 8, 9, 16, 17, 9, 10, 1, 1, -2, 2, 3, 3, 14, 4, 5, 10, 13, 9, 9, 10, 1, 0, 0, 5, 5, 6, 6, 13, 13, 14, -14, 4, 4, 5, 11, 11, 12, 15, 15, 16}; 1642550becSJunchao Zhang PetscInt j0[] = {1, 6, 2, 4, 10, 15, 13, 16, 11, 2, 7, 3, 8, 4, 9, 13, 5, 2, 15, 14, 10, 16, 11, 2, 0, 1, 6, -11, 0, 7, 15, 17, 11, 13, 5, 8, 2, 12, 17, 13, 3, 16, 9}; 1742550becSJunchao Zhang 1842550becSJunchao Zhang PetscInt i1[] = {8, 5, 15, 16, 6, 13, 4, 17, 8, 9, 9, 10, -6, 12, 7, 3, -4, 1, 1, 2, 5, 5, 6, 14, 17, 8, 9, 9, 10, 4, 5, 10, 11, 1, 2}; 1942550becSJunchao Zhang PetscInt j1[] = {2, 3, 16, 9, 5, 17, 1, 13, 4, 10, 16, 11, -5, 12, 1, 7, -1, 2, 7, 3, 6, 11, 0, 11, 13, 4, 10, 16, 11, 8, -2, 15, 12, 7, 3}; 2042550becSJunchao Zhang 2142550becSJunchao Zhang PetscInt i2[] = {3, 4, 1, 10, 0, 1, 1, 2, 1, 1, 2, 2, 3, 3, 4, 4, 1, 2, 5, 5, 6, 4, 17, 0, 1, 1, 8, 5, 5, 6, 4, 7, 8, 5}; 2242550becSJunchao Zhang PetscInt j2[] = {7, 1, 2, 11, 5, 2, 7, 3, 2, 7, 3, 8, 4, 9, 3, 5, 7, 3, 6, 11, 0, 1, 13, 5, 2, 7, 4, 6, 11, 0, 1, 3, 4, 2}; 2342550becSJunchao Zhang 2442550becSJunchao Zhang struct { 2542550becSJunchao Zhang PetscInt *i,*j,n; 2642550becSJunchao Zhang } coo[3] = {{i0,j0,sizeof(i0)/sizeof(PetscInt)}, {i1,j1,sizeof(i1)/sizeof(PetscInt)}, {i2,j2,sizeof(i2)/sizeof(PetscInt)}}; 2742550becSJunchao Zhang 28*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-ignore_remote",&flg,NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 315f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 3242550becSJunchao Zhang 332c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 3,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"This test requires at most 3 processes"); 3442550becSJunchao Zhang 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 36*b122ec5aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 38*b122ec5aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,2,NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,2,NULL,2,NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_IGNORE_OFF_PROC_ENTRIES,flg)); 4242550becSJunchao Zhang 4342550becSJunchao Zhang for (k=0; k<coo[rank].n; k++) { 4442550becSJunchao Zhang PetscScalar val = coo[rank].j[k]; 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A,coo[rank].i[k],coo[rank].j[k],val,ADD_VALUES)); 4642550becSJunchao Zhang } 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 4942550becSJunchao Zhang 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 51*b122ec5aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_IGNORE_OFF_PROC_ENTRIES,flg)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetPreallocationCOO(B,coo[rank].n,coo[rank].i,coo[rank].j)); 5542550becSJunchao Zhang 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(coo[rank].n,&vals)); 5742550becSJunchao Zhang for (k=0; k<coo[rank].n; k++) vals[k] = coo[rank].j[k]; 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesCOO(B,vals,ADD_VALUES)); 5942550becSJunchao Zhang 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(A,B,&equal)); 6142550becSJunchao Zhang 6242550becSJunchao Zhang if (!equal) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 6542550becSJunchao Zhang SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"MatSetValuesCOO() failed"); 6642550becSJunchao Zhang } 6742550becSJunchao Zhang 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vals)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 7142550becSJunchao Zhang 72*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 73*b122ec5aSJacob Faibussowitsch return 0; 7442550becSJunchao Zhang } 7542550becSJunchao Zhang 7642550becSJunchao Zhang /*TEST 7742550becSJunchao Zhang 7821c60b78SJunchao Zhang testset: 7921c60b78SJunchao Zhang output_file: output/ex254_1.out 8042550becSJunchao Zhang nsize: {{1 2 3}} 8121c60b78SJunchao Zhang args: -ignore_remote {{0 1}} 8221c60b78SJunchao Zhang 8321c60b78SJunchao Zhang test: 8421c60b78SJunchao Zhang suffix: kokkos 8542550becSJunchao Zhang requires: kokkos_kernels 8642550becSJunchao Zhang args: -mat_type aijkokkos 8721c60b78SJunchao Zhang 8821c60b78SJunchao Zhang test: 8921c60b78SJunchao Zhang suffix: cuda 9021c60b78SJunchao Zhang requires: cuda 9121c60b78SJunchao Zhang args: -mat_type aijcusparse 9221c60b78SJunchao Zhang 9321c60b78SJunchao Zhang test: 9421c60b78SJunchao Zhang suffix: aij 9521c60b78SJunchao Zhang args: -mat_type aij 9642550becSJunchao Zhang 9742550becSJunchao Zhang TEST*/ 98