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