xref: /petsc/src/mat/tests/ex260.c (revision 38b83642f90ca7c29d17300cee5f7fcf84bc0cc1)
1*38b83642SBarry Smith static char help[] = "Tests that MatView() and MatLoad() work for MPIAIJ matrix with total nz > PETSC_MAX_INT\n\n";
2*38b83642SBarry Smith 
3*38b83642SBarry Smith #include <petscmat.h>
4*38b83642SBarry Smith 
5*38b83642SBarry Smith int main(int argc, char **args)
6*38b83642SBarry Smith {
7*38b83642SBarry Smith   PetscInt     n = 1000000, nzr = (PetscInt)((double)PETSC_MAX_INT) / (3.8 * n);
8*38b83642SBarry Smith   Mat          A;
9*38b83642SBarry Smith   PetscScalar *a;
10*38b83642SBarry Smith   PetscInt    *ii, *jd, *jo;
11*38b83642SBarry Smith   PetscMPIInt  rank, size;
12*38b83642SBarry Smith   PetscViewer  viewer;
13*38b83642SBarry Smith 
14*38b83642SBarry Smith   PetscFunctionBeginUser;
15*38b83642SBarry Smith   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
16*38b83642SBarry Smith 
17*38b83642SBarry Smith   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18*38b83642SBarry Smith   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with two MPI ranks");
19*38b83642SBarry Smith   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20*38b83642SBarry Smith   PetscCall(PetscMalloc4(n + 1, &ii, n * nzr, &jd, n * nzr, &jo, n * nzr, &a));
21*38b83642SBarry Smith   ii[0] = 0;
22*38b83642SBarry Smith   for (PetscInt i = 0; i < n; i++) {
23*38b83642SBarry Smith     ii[i + 1] = ii[i] + nzr;
24*38b83642SBarry Smith     for (PetscInt j = 0; j < nzr; j++) jd[ii[i] + j] = j;
25*38b83642SBarry Smith     if (rank == 0) {
26*38b83642SBarry Smith       for (PetscInt j = 0; j < nzr; j++) jo[ii[i] + j] = n + j - 1;
27*38b83642SBarry Smith     } else {
28*38b83642SBarry Smith       for (PetscInt j = 0; j < nzr; j++) jo[ii[i] + j] = j;
29*38b83642SBarry Smith     }
30*38b83642SBarry Smith   }
31*38b83642SBarry Smith   PetscCall(MatCreateMPIAIJWithSplitArrays(PETSC_COMM_WORLD, n, n, PETSC_DETERMINE, PETSC_DETERMINE, ii, jd, a, ii, jo, a, &A));
32*38b83642SBarry Smith   PetscCall(MatView(A, PETSC_VIEWER_BINARY_WORLD));
33*38b83642SBarry Smith   PetscCall(MatDestroy(&A));
34*38b83642SBarry Smith   PetscCall(PetscFree4(ii, jd, jo, a));
35*38b83642SBarry Smith 
36*38b83642SBarry Smith   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "binaryoutput", FILE_MODE_READ, &viewer));
37*38b83642SBarry Smith   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
38*38b83642SBarry Smith   PetscCall(MatLoad(A, viewer));
39*38b83642SBarry Smith   PetscCall(PetscViewerDestroy(&viewer));
40*38b83642SBarry Smith   PetscCall(MatDestroy(&A));
41*38b83642SBarry Smith   PetscCall(PetscFinalize());
42*38b83642SBarry Smith   return 0;
43*38b83642SBarry Smith }
44*38b83642SBarry Smith 
45*38b83642SBarry Smith /*TEST
46*38b83642SBarry Smith 
47*38b83642SBarry Smith    test:
48*38b83642SBarry Smith       TODO: requires to much memory to run in the CI
49*38b83642SBarry Smith       nsize: 2
50*38b83642SBarry Smith 
51*38b83642SBarry Smith TEST*/
52