Actual source code: aijhdf5.c
1: #include <../src/mat/impls/aij/seq/aij.h>
2: #include <petsc/private/isimpl.h>
3: #include <petsc/private/vecimpl.h>
4: #include <petsclayouthdf5.h>
6: #if defined(PETSC_HAVE_HDF5)
7: PetscErrorCode MatLoad_AIJ_HDF5(Mat mat, PetscViewer viewer)
8: {
9: PetscViewerFormat format;
10: const PetscInt *i_glob = NULL;
11: PetscInt *i = NULL;
12: const PetscInt *j = NULL;
13: const PetscScalar *a = NULL;
14: char *a_name = NULL, *i_name = NULL, *j_name = NULL, *c_name = NULL;
15: const char *mat_name = NULL;
16: PetscInt p, m, M, N;
17: PetscInt bs = mat->rmap->bs;
18: PetscInt *range;
19: PetscBool flg;
20: IS is_i = NULL, is_j = NULL;
21: Vec vec_a = NULL;
22: PetscLayout jmap = NULL;
23: MPI_Comm comm;
24: PetscMPIInt rank, size;
26: PetscFunctionBegin;
27: PetscCall(PetscViewerGetFormat(viewer, &format));
28: switch (format) {
29: case PETSC_VIEWER_HDF5_PETSC:
30: case PETSC_VIEWER_DEFAULT:
31: case PETSC_VIEWER_NATIVE:
32: case PETSC_VIEWER_HDF5_MAT:
33: break;
34: default:
35: SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
36: }
38: PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
39: PetscCallMPI(MPI_Comm_rank(comm, &rank));
40: PetscCallMPI(MPI_Comm_size(comm, &size));
41: PetscCall(PetscObjectGetName((PetscObject)mat, &mat_name));
42: if (format == PETSC_VIEWER_HDF5_MAT) {
43: PetscCall(PetscStrallocpy("jc", &i_name));
44: PetscCall(PetscStrallocpy("ir", &j_name));
45: PetscCall(PetscStrallocpy("data", &a_name));
46: PetscCall(PetscStrallocpy("MATLAB_sparse", &c_name));
47: } else {
48: /* TODO Once corresponding MatView is implemented, change the names to i,j,a */
49: /* TODO Maybe there could be both namings in the file, using "symbolic link" features of HDF5. */
50: PetscCall(PetscStrallocpy("jc", &i_name));
51: PetscCall(PetscStrallocpy("ir", &j_name));
52: PetscCall(PetscStrallocpy("data", &a_name));
53: PetscCall(PetscStrallocpy("MATLAB_sparse", &c_name));
54: }
56: PetscOptionsBegin(comm, NULL, "Options for loading matrix from HDF5", "Mat");
57: PetscCall(PetscOptionsInt("-matload_block_size", "Set the blocksize used to store the matrix", "MatLoad", bs, &bs, &flg));
58: PetscOptionsEnd();
59: if (flg) PetscCall(MatSetBlockSize(mat, bs));
61: PetscCall(PetscViewerHDF5PushGroup(viewer, mat_name));
62: PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, c_name, PETSC_INT, NULL, &N));
63: PetscCall(PetscViewerHDF5ReadSizes(viewer, i_name, NULL, &M));
64: --M; /* i has size M+1 as there is global number of nonzeros stored at the end */
66: if (format == PETSC_VIEWER_HDF5_MAT && mat->symmetric != PETSC_BOOL3_TRUE) {
67: /* Swap row and columns layout for unallocated matrix. I want to avoid calling MatTranspose() just to transpose sparsity pattern and layout. */
68: if (!mat->preallocated) {
69: PetscLayout tmp;
70: tmp = mat->rmap;
71: mat->rmap = mat->cmap;
72: mat->cmap = tmp;
73: } else SETERRQ(comm, PETSC_ERR_SUP, "Not for preallocated matrix - we would need to transpose it here which we want to avoid");
74: }
76: /* If global sizes are set, check if they are consistent with that given in the file */
77: PetscCheck(mat->rmap->N < 0 || mat->rmap->N == M, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows: Matrix in file has (%" PetscInt_FMT ") and input matrix has (%" PetscInt_FMT ")", mat->rmap->N, M);
78: PetscCheck(mat->cmap->N < 0 || mat->cmap->N == N, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols: Matrix in file has (%" PetscInt_FMT ") and input matrix has (%" PetscInt_FMT ")", mat->cmap->N, N);
80: /* Determine ownership of all (block) rows and columns */
81: mat->rmap->N = M;
82: mat->cmap->N = N;
83: PetscCall(PetscLayoutSetUp(mat->rmap));
84: PetscCall(PetscLayoutSetUp(mat->cmap));
85: m = mat->rmap->n;
87: /* Read array i (array of row indices) */
88: PetscCall(PetscMalloc1(m + 1, &i)); /* allocate i with one more position for local number of nonzeros on each rank */
89: i[0] = i[m] = 0; /* make the last entry always defined - the code block below overwrites it just on last rank */
90: if (rank == size - 1) m++; /* in the loaded array i_glob, only the last rank has one more position with the global number of nonzeros */
91: M++;
92: PetscCall(ISCreate(comm, &is_i));
93: PetscCall(PetscObjectSetName((PetscObject)is_i, i_name));
94: PetscCall(PetscLayoutSetLocalSize(is_i->map, m));
95: PetscCall(PetscLayoutSetSize(is_i->map, M));
96: PetscCall(ISLoad(is_i, viewer));
97: PetscCall(ISGetIndices(is_i, &i_glob));
98: PetscCall(PetscArraycpy(i, i_glob, m));
100: /* Reset m and M to the matrix sizes */
101: m = mat->rmap->n;
102: M--;
104: /* Create PetscLayout for j and a vectors; construct ranges first */
105: PetscCall(PetscMalloc1(size + 1, &range));
106: PetscCallMPI(MPI_Allgather(i, 1, MPIU_INT, range, 1, MPIU_INT, comm));
107: /* Last rank has global number of nonzeros (= length of j and a arrays) in i[m] (last i entry) so broadcast it */
108: range[size] = i[m];
109: PetscCallMPI(MPI_Bcast(&range[size], 1, MPIU_INT, size - 1, comm));
110: for (p = size - 1; p > 0; p--) {
111: if (!range[p]) range[p] = range[p + 1]; /* for ranks with 0 rows, take the value from the next processor */
112: }
113: i[m] = range[rank + 1]; /* i[m] (last i entry) is equal to next rank's offset */
114: /* Deduce rstart, rend, n and N from the ranges */
115: PetscCall(PetscLayoutCreateFromRanges(comm, range, PETSC_OWN_POINTER, 1, &jmap));
117: /* Convert global to local indexing of rows */
118: for (p = 1; p < m + 1; ++p) i[p] -= i[0];
119: i[0] = 0;
121: /* Read array j (array of column indices) */
122: PetscCall(ISCreate(comm, &is_j));
123: PetscCall(PetscObjectSetName((PetscObject)is_j, j_name));
124: PetscCall(PetscLayoutDuplicate(jmap, &is_j->map));
125: PetscCall(ISLoad(is_j, viewer));
126: PetscCall(ISGetIndices(is_j, &j));
128: /* Read array a (array of values) */
129: PetscCall(VecCreate(comm, &vec_a));
130: PetscCall(PetscObjectSetName((PetscObject)vec_a, a_name));
131: PetscCall(PetscLayoutDuplicate(jmap, &vec_a->map));
132: PetscCall(VecLoad(vec_a, viewer));
133: PetscCall(VecGetArrayRead(vec_a, &a));
135: /* populate matrix */
136: if (!((PetscObject)mat)->type_name) PetscCall(MatSetType(mat, MATAIJ));
137: PetscCall(MatSeqAIJSetPreallocationCSR(mat, i, j, a));
138: PetscCall(MatMPIAIJSetPreallocationCSR(mat, i, j, a));
139: /*
140: PetscCall(MatSeqBAIJSetPreallocationCSR(mat,bs,i,j,a));
141: PetscCall(MatMPIBAIJSetPreallocationCSR(mat,bs,i,j,a));
142: */
144: if (format == PETSC_VIEWER_HDF5_MAT && mat->symmetric != PETSC_BOOL3_TRUE) {
145: /* Transpose the input matrix back */
146: PetscCall(MatTranspose(mat, MAT_INPLACE_MATRIX, &mat));
147: }
149: PetscCall(PetscViewerHDF5PopGroup(viewer));
150: PetscCall(PetscFree(i_name));
151: PetscCall(PetscFree(j_name));
152: PetscCall(PetscFree(a_name));
153: PetscCall(PetscFree(c_name));
154: PetscCall(PetscLayoutDestroy(&jmap));
155: PetscCall(PetscFree(i));
156: PetscCall(ISRestoreIndices(is_i, &i_glob));
157: PetscCall(ISRestoreIndices(is_j, &j));
158: PetscCall(VecRestoreArrayRead(vec_a, &a));
159: PetscCall(ISDestroy(&is_i));
160: PetscCall(ISDestroy(&is_j));
161: PetscCall(VecDestroy(&vec_a));
162: PetscFunctionReturn(PETSC_SUCCESS);
163: }
164: #endif