xref: /petsc/src/mat/tests/ex254.c (revision 42550bec5353b99bbf937a3bdb81525a3b04a032)
1*42550becSJunchao Zhang static char help[] = "Test MatSetValuesCOO for MPIAIJKOKKOS mat \n\n";
2*42550becSJunchao Zhang 
3*42550becSJunchao Zhang #include <petscmat.h>
4*42550becSJunchao Zhang int main(int argc,char **args)
5*42550becSJunchao Zhang {
6*42550becSJunchao Zhang   PetscErrorCode  ierr;
7*42550becSJunchao Zhang   Mat             A,B;
8*42550becSJunchao Zhang   PetscInt        k;
9*42550becSJunchao Zhang   const PetscInt  M = 18,N = 18;
10*42550becSJunchao Zhang   PetscMPIInt     rank,size;
11*42550becSJunchao Zhang   PetscBool       equal;
12*42550becSJunchao Zhang   PetscScalar     *vals;
13*42550becSJunchao Zhang 
14*42550becSJunchao Zhang   /* Construct 18 x 18 matrices, which are big enough to have complex communication patterns but still small enough for debugging */
15*42550becSJunchao 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};
16*42550becSJunchao 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};
17*42550becSJunchao Zhang 
18*42550becSJunchao 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};
19*42550becSJunchao 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};
20*42550becSJunchao Zhang 
21*42550becSJunchao 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};
22*42550becSJunchao 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};
23*42550becSJunchao Zhang 
24*42550becSJunchao Zhang   struct {
25*42550becSJunchao Zhang     PetscInt *i,*j,n;
26*42550becSJunchao Zhang   } coo[3] = {{i0,j0,sizeof(i0)/sizeof(PetscInt)}, {i1,j1,sizeof(i1)/sizeof(PetscInt)}, {i2,j2,sizeof(i2)/sizeof(PetscInt)}};
27*42550becSJunchao Zhang 
28*42550becSJunchao Zhang   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29*42550becSJunchao Zhang   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
30*42550becSJunchao Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
31*42550becSJunchao Zhang 
32*42550becSJunchao Zhang   if (size > 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"This test requires at most 3 processes");
33*42550becSJunchao Zhang 
34*42550becSJunchao Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
35*42550becSJunchao Zhang   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);
36*42550becSJunchao Zhang   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
37*42550becSJunchao Zhang   ierr = MatSeqAIJSetPreallocation(A,2,NULL);
38*42550becSJunchao Zhang   ierr = MatMPIAIJSetPreallocation(A,2,NULL,2,NULL);CHKERRQ(ierr);
39*42550becSJunchao Zhang   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
40*42550becSJunchao Zhang 
41*42550becSJunchao Zhang   for (k=0; k<coo[rank].n; k++) {
42*42550becSJunchao Zhang     PetscScalar val = coo[rank].j[k];
43*42550becSJunchao Zhang     ierr = MatSetValue(A,coo[rank].i[k],coo[rank].j[k],val,ADD_VALUES);CHKERRQ(ierr);
44*42550becSJunchao Zhang   }
45*42550becSJunchao Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46*42550becSJunchao Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47*42550becSJunchao Zhang 
48*42550becSJunchao Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
49*42550becSJunchao Zhang   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
50*42550becSJunchao Zhang   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
51*42550becSJunchao Zhang   ierr = MatSetPreallocationCOO(B,coo[rank].n,coo[rank].i,coo[rank].j);CHKERRQ(ierr);
52*42550becSJunchao Zhang 
53*42550becSJunchao Zhang   ierr = PetscMalloc1(coo[rank].n,&vals);CHKERRQ(ierr);
54*42550becSJunchao Zhang   for (k=0; k<coo[rank].n; k++) vals[k] = coo[rank].j[k];
55*42550becSJunchao Zhang   ierr = MatSetValuesCOO(B,vals,ADD_VALUES);CHKERRQ(ierr);
56*42550becSJunchao Zhang 
57*42550becSJunchao Zhang   ierr = MatEqual(A,B,&equal);CHKERRQ(ierr);
58*42550becSJunchao Zhang 
59*42550becSJunchao Zhang   if (!equal) {
60*42550becSJunchao Zhang     ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
61*42550becSJunchao Zhang     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
62*42550becSJunchao Zhang     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"MatSetValuesCOO() failed");
63*42550becSJunchao Zhang   }
64*42550becSJunchao Zhang 
65*42550becSJunchao Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
66*42550becSJunchao Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
67*42550becSJunchao Zhang   ierr = MatDestroy(&B);CHKERRQ(ierr);
68*42550becSJunchao Zhang 
69*42550becSJunchao Zhang   ierr = PetscFinalize();
70*42550becSJunchao Zhang   return ierr;
71*42550becSJunchao Zhang }
72*42550becSJunchao Zhang 
73*42550becSJunchao Zhang /*TEST
74*42550becSJunchao Zhang 
75*42550becSJunchao Zhang   test:
76*42550becSJunchao Zhang     nsize: {{1 2 3}}
77*42550becSJunchao Zhang     requires: kokkos_kernels
78*42550becSJunchao Zhang     args: -mat_type aijkokkos
79*42550becSJunchao Zhang     output_file: output/ex254_1.out
80*42550becSJunchao Zhang 
81*42550becSJunchao Zhang TEST*/
82*42550becSJunchao Zhang 
83