xref: /petsc/src/mat/tests/ex254.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
142550becSJunchao Zhang static char help[] = "Test MatSetValuesCOO for MPIAIJKOKKOS mat \n\n";
242550becSJunchao Zhang 
342550becSJunchao Zhang #include <petscmat.h>
442550becSJunchao Zhang int main(int argc,char **args)
542550becSJunchao Zhang {
642550becSJunchao Zhang   PetscErrorCode  ierr;
742550becSJunchao Zhang   Mat             A,B;
842550becSJunchao Zhang   PetscInt        k;
942550becSJunchao Zhang   const PetscInt  M = 18,N = 18;
1042550becSJunchao Zhang   PetscMPIInt     rank,size;
1142550becSJunchao Zhang   PetscBool       equal;
1242550becSJunchao Zhang   PetscScalar     *vals;
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 
2842550becSJunchao Zhang   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
2942550becSJunchao Zhang   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
3042550becSJunchao Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
3142550becSJunchao Zhang 
32*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 3,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"This test requires at most 3 processes");
3342550becSJunchao Zhang 
3442550becSJunchao Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
3542550becSJunchao Zhang   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);
3642550becSJunchao Zhang   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
3742550becSJunchao Zhang   ierr = MatSeqAIJSetPreallocation(A,2,NULL);
3842550becSJunchao Zhang   ierr = MatMPIAIJSetPreallocation(A,2,NULL,2,NULL);CHKERRQ(ierr);
3942550becSJunchao Zhang   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
4042550becSJunchao Zhang 
4142550becSJunchao Zhang   for (k=0; k<coo[rank].n; k++) {
4242550becSJunchao Zhang     PetscScalar val = coo[rank].j[k];
4342550becSJunchao Zhang     ierr = MatSetValue(A,coo[rank].i[k],coo[rank].j[k],val,ADD_VALUES);CHKERRQ(ierr);
4442550becSJunchao Zhang   }
4542550becSJunchao Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4642550becSJunchao Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4742550becSJunchao Zhang 
4842550becSJunchao Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
4942550becSJunchao Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
5042550becSJunchao Zhang   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
5142550becSJunchao Zhang   ierr = MatSetPreallocationCOO(B,coo[rank].n,coo[rank].i,coo[rank].j);CHKERRQ(ierr);
5242550becSJunchao Zhang 
5342550becSJunchao Zhang   ierr = PetscMalloc1(coo[rank].n,&vals);CHKERRQ(ierr);
5442550becSJunchao Zhang   for (k=0; k<coo[rank].n; k++) vals[k] = coo[rank].j[k];
5542550becSJunchao Zhang   ierr = MatSetValuesCOO(B,vals,ADD_VALUES);CHKERRQ(ierr);
5642550becSJunchao Zhang 
5742550becSJunchao Zhang   ierr = MatEqual(A,B,&equal);CHKERRQ(ierr);
5842550becSJunchao Zhang 
5942550becSJunchao Zhang   if (!equal) {
6042550becSJunchao Zhang     ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
6142550becSJunchao Zhang     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
6242550becSJunchao Zhang     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"MatSetValuesCOO() failed");
6342550becSJunchao Zhang   }
6442550becSJunchao Zhang 
6542550becSJunchao Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
6642550becSJunchao Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
6742550becSJunchao Zhang   ierr = MatDestroy(&B);CHKERRQ(ierr);
6842550becSJunchao Zhang 
6942550becSJunchao Zhang   ierr = PetscFinalize();
7042550becSJunchao Zhang   return ierr;
7142550becSJunchao Zhang }
7242550becSJunchao Zhang 
7342550becSJunchao Zhang /*TEST
7442550becSJunchao Zhang 
7542550becSJunchao Zhang   test:
7642550becSJunchao Zhang     nsize: {{1 2 3}}
7742550becSJunchao Zhang     requires: kokkos_kernels
7842550becSJunchao Zhang     args: -mat_type aijkokkos
7942550becSJunchao Zhang     output_file: output/ex254_1.out
8042550becSJunchao Zhang 
8142550becSJunchao Zhang TEST*/
8242550becSJunchao Zhang 
83