xref: /petsc/src/mat/tests/ex303.c (revision f0c227fb9e7ece2f21a83203a567dc40969a87c2)
1c0c276a7Ssdargavi static char help[] = "Testing MatCreateMPIAIJWithSeqAIJ().\n\n";
2c0c276a7Ssdargavi 
3c0c276a7Ssdargavi #include <petscmat.h>
4c0c276a7Ssdargavi 
5c0c276a7Ssdargavi int main(int argc, char **argv)
6c0c276a7Ssdargavi {
7*f0c227fbSJunchao Zhang   Mat                mat, mat2;
8*f0c227fbSJunchao Zhang   Mat                A, B;   // diag, offdiag of mat
9*f0c227fbSJunchao Zhang   Mat                A2, B2; // diag, offdiag of mat2
10*f0c227fbSJunchao Zhang   PetscInt           M, N, m = 5, n = 6, d_nz = 3, o_nz = 4;
11*f0c227fbSJunchao Zhang   PetscInt          *bi, *bj, nd;
12*f0c227fbSJunchao Zhang   const PetscScalar *ba;
13c0c276a7Ssdargavi   const PetscInt    *garray;
14c0c276a7Ssdargavi   PetscInt          *garray_h;
15c0c276a7Ssdargavi   PetscBool          equal, done;
16*f0c227fbSJunchao Zhang   PetscMPIInt        size;
17*f0c227fbSJunchao Zhang   char               mat_type[PETSC_MAX_PATH_LEN];
18*f0c227fbSJunchao Zhang   PetscBool          flg;
19c0c276a7Ssdargavi 
20c0c276a7Ssdargavi   PetscFunctionBeginUser;
21c0c276a7Ssdargavi   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
22c0c276a7Ssdargavi   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23c0c276a7Ssdargavi   PetscCheck(size > 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 2 or more processes");
24c0c276a7Ssdargavi 
25*f0c227fbSJunchao Zhang   // Create a MATMPIAIJ matrix for checking against
26*f0c227fbSJunchao Zhang   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, m, n, PETSC_DECIDE, PETSC_DECIDE, d_nz, NULL, o_nz, NULL, &mat));
27*f0c227fbSJunchao Zhang   PetscCall(MatSetRandom(mat, NULL));
28*f0c227fbSJunchao Zhang   PetscCall(MatGetSize(mat, &M, &N));
29*f0c227fbSJunchao Zhang   PetscCall(MatGetLocalSize(mat, &m, &n));
30*f0c227fbSJunchao Zhang   PetscCall(MatMPIAIJGetSeqAIJ(mat, &A, &B, &garray));
31c0c276a7Ssdargavi 
32*f0c227fbSJunchao Zhang   PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_type", mat_type, sizeof(mat_type), &flg));
33*f0c227fbSJunchao Zhang   if (!flg) PetscCall(PetscStrncpy(mat_type, MATSEQAIJ, sizeof(mat_type))); // Default to MATAIJ
34*f0c227fbSJunchao Zhang   PetscCall(MatConvert(A, mat_type, MAT_INITIAL_MATRIX, &A2));              // Copy A, B to A2, B2 but in the given mat_type
35*f0c227fbSJunchao Zhang   PetscCall(MatConvert(B, mat_type, MAT_INITIAL_MATRIX, &B2));
36c0c276a7Ssdargavi 
37*f0c227fbSJunchao Zhang   // The garray passed in has to be on the host, but it can be created on device and copied to the host.
38*f0c227fbSJunchao Zhang   // We're just going to copy the existing host values here.
39c0c276a7Ssdargavi   PetscInt nonlocalCols;
40*f0c227fbSJunchao Zhang   PetscCall(MatGetLocalSize(B, NULL, &nonlocalCols));
41c0c276a7Ssdargavi   PetscCall(PetscMalloc1(nonlocalCols, &garray_h));
42c0c276a7Ssdargavi   for (int i = 0; i < nonlocalCols; i++) { garray_h[i] = garray[i]; }
43c0c276a7Ssdargavi 
44*f0c227fbSJunchao Zhang   // 1) Test MatCreateMPIAIJWithSeqAIJ with compressed off-diag.
45*f0c227fbSJunchao Zhang   PetscCall(MatCreateMPIAIJWithSeqAIJ(PETSC_COMM_WORLD, M, N, A2, B2, garray_h, &mat2));
46*f0c227fbSJunchao Zhang   PetscCall(MatEqual(mat, mat2, &equal));                   // might directly compare value arrays
47*f0c227fbSJunchao Zhang   if (equal) PetscCall(MatMultEqual(mat, mat2, 2, &equal)); // compare via MatMult()
48c0c276a7Ssdargavi   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Likely a bug in MatCreateMPIAIJWithSeqAIJ()");
49*f0c227fbSJunchao Zhang   PetscCall(MatDestroy(&mat2)); // A2 and B2 are also destroyed
50c0c276a7Ssdargavi 
51*f0c227fbSJunchao Zhang   // 2) Test MatCreateMPIAIJWithSeqAIJ with non-compressed off-diag.
52*f0c227fbSJunchao Zhang   PetscCall(MatConvert(A, mat_type, MAT_INITIAL_MATRIX, &A2));
53c0c276a7Ssdargavi 
54*f0c227fbSJunchao Zhang   // Create a version of B2 of size N with global indices
55*f0c227fbSJunchao Zhang   PetscCall(MatGetRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&bi, (const PetscInt **)&bj, &done));
56*f0c227fbSJunchao Zhang   PetscCall(MatSeqAIJGetArrayRead(B, &ba));
57*f0c227fbSJunchao Zhang   PetscCall(MatCreate(PETSC_COMM_SELF, &B2));
58*f0c227fbSJunchao Zhang   PetscCall(MatSetSizes(B2, m, N, m, N));
59*f0c227fbSJunchao Zhang   PetscCall(MatSetType(B2, mat_type));
60*f0c227fbSJunchao Zhang   PetscCall(MatSeqAIJSetPreallocation(B2, o_nz, NULL));
61*f0c227fbSJunchao Zhang   for (int i = 0; i < m; i++) { // Fill B2 with values from B
62*f0c227fbSJunchao Zhang     for (int j = 0; j < bi[i + 1] - bi[i]; j++) { PetscCall(MatSetValue(B2, i, garray[bj[bi[i] + j]], ba[bi[i] + j], INSERT_VALUES)); }
63c0c276a7Ssdargavi   }
64*f0c227fbSJunchao Zhang   PetscCall(MatAssemblyBegin(B2, MAT_FINAL_ASSEMBLY));
65*f0c227fbSJunchao Zhang   PetscCall(MatAssemblyEnd(B2, MAT_FINAL_ASSEMBLY));
66c0c276a7Ssdargavi 
67*f0c227fbSJunchao Zhang   PetscCall(MatRestoreRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nd, (const PetscInt **)&bi, (const PetscInt **)&bj, &done));
68*f0c227fbSJunchao Zhang   PetscCall(MatSeqAIJRestoreArrayRead(B, &ba));
69c0c276a7Ssdargavi 
70*f0c227fbSJunchao Zhang   PetscCall(MatCreateMPIAIJWithSeqAIJ(PETSC_COMM_WORLD, M, N, A2, B2, NULL, &mat2));
71*f0c227fbSJunchao Zhang   PetscCall(MatEqual(mat, mat2, &equal));
72*f0c227fbSJunchao Zhang   if (equal) PetscCall(MatMultEqual(mat, mat2, 2, &equal));
73c0c276a7Ssdargavi   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Likely a bug in MatCreateMPIAIJWithSeqAIJ()");
74c0c276a7Ssdargavi 
75*f0c227fbSJunchao Zhang   PetscCall(MatDestroy(&mat));
76*f0c227fbSJunchao Zhang   PetscCall(MatDestroy(&mat2));
77c0c276a7Ssdargavi   PetscCall(PetscFinalize());
78c0c276a7Ssdargavi   return 0;
79c0c276a7Ssdargavi }
80c0c276a7Ssdargavi 
81c0c276a7Ssdargavi /*TEST
82c0c276a7Ssdargavi 
83*f0c227fbSJunchao Zhang   testset:
84c0c276a7Ssdargavi     nsize: 2
85c0c276a7Ssdargavi     output_file: output/empty.out
86c0c276a7Ssdargavi 
87*f0c227fbSJunchao Zhang     test:
88*f0c227fbSJunchao Zhang       suffix: 1
89*f0c227fbSJunchao Zhang       args: -mat_type aij
90*f0c227fbSJunchao Zhang 
91*f0c227fbSJunchao Zhang     test:
92*f0c227fbSJunchao Zhang       suffix: 2
93*f0c227fbSJunchao Zhang       args: -mat_type aijkokkos
94*f0c227fbSJunchao Zhang       requires: kokkos_kernels
95*f0c227fbSJunchao Zhang 
96*f0c227fbSJunchao Zhang     test:
97*f0c227fbSJunchao Zhang       suffix: 3
98*f0c227fbSJunchao Zhang       args: -mat_type aijcusparse
99*f0c227fbSJunchao Zhang       requires: cuda
100*f0c227fbSJunchao Zhang 
101*f0c227fbSJunchao Zhang     test:
102*f0c227fbSJunchao Zhang       suffix: 4
103*f0c227fbSJunchao Zhang       args: -mat_type aijhipsparse
104*f0c227fbSJunchao Zhang       requires: hip
105*f0c227fbSJunchao Zhang 
106c0c276a7Ssdargavi TEST*/
107