Actual source code: aijhdf5.c

petsc-3.12.5 2020-03-29
Report Typos and Errors

  2:  #include <../src/mat/impls/aij/seq/aij.h>
  3:  #include <petsc/private/isimpl.h>
  4:  #include <petsc/private/vecimpl.h>
  5: #include <petsclayouthdf5.h>

  7: #if defined(PETSC_HAVE_HDF5)
  8: PetscErrorCode MatLoad_AIJ_HDF5(Mat mat, PetscViewer viewer)
  9: {
 10:   PetscViewerFormat format;
 11:   const PetscInt  *i_glob = NULL;
 12:   PetscInt        *i = NULL;
 13:   const PetscInt  *j = NULL;
 14:   const PetscScalar *a = NULL;
 15:   char            *a_name = NULL, *i_name = NULL, *j_name = NULL, *c_name = NULL;
 16:   const char      *mat_name = NULL;
 17:   PetscInt        p, m, M, N;
 18:   PetscInt        bs = mat->rmap->bs;
 19:   PetscInt        *range;
 20:   PetscBool       flg;
 21:   IS              is_i = NULL, is_j = NULL;
 22:   Vec             vec_a = NULL;
 23:   PetscLayout     jmap = NULL;
 24:   MPI_Comm        comm;
 25:   PetscMPIInt     rank, size;
 26:   PetscErrorCode  ierr;

 29:   PetscViewerGetFormat(viewer, &format);
 30:   switch (format) {
 31:     case PETSC_VIEWER_HDF5_PETSC:
 32:     case PETSC_VIEWER_DEFAULT:
 33:     case PETSC_VIEWER_NATIVE:
 34:     case PETSC_VIEWER_HDF5_MAT:
 35:       break;
 36:     default:
 37:       SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"PetscViewerFormat %s not supported for HDF5 input.",PetscViewerFormats[format]);
 38:   }

 40:   PetscObjectGetComm((PetscObject)mat,&comm);
 41:   MPI_Comm_rank(comm,&rank);
 42:   MPI_Comm_size(comm,&size);
 43:   PetscObjectGetName((PetscObject)mat,&mat_name);
 44:   if (format==PETSC_VIEWER_HDF5_MAT) {
 45:     PetscStrallocpy("jc",&i_name);
 46:     PetscStrallocpy("ir",&j_name);
 47:     PetscStrallocpy("data",&a_name);
 48:     PetscStrallocpy("MATLAB_sparse",&c_name);
 49:   } else {
 50:     /* TODO Once corresponding MatView is implemented, change the names to i,j,a */
 51:     /* TODO Maybe there could be both namings in the file, using "symbolic link" features of HDF5. */
 52:     PetscStrallocpy("jc",&i_name);
 53:     PetscStrallocpy("ir",&j_name);
 54:     PetscStrallocpy("data",&a_name);
 55:     PetscStrallocpy("MATLAB_sparse",&c_name);
 56:   }

 58:   PetscOptionsBegin(comm,NULL,"Options for loading matrix from HDF5","Mat");
 59:   PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,&flg);
 60:   PetscOptionsEnd();
 61:   if (flg) {
 62:     MatSetBlockSize(mat, bs);
 63:   }

 65:   PetscViewerHDF5PushGroup(viewer,mat_name);
 66:   PetscViewerHDF5ReadAttribute(viewer,NULL,c_name,PETSC_INT,&N);
 67:   PetscViewerHDF5ReadSizes(viewer, i_name, NULL, &M);
 68:   --M;  /* i has size M+1 as there is global number of nonzeros stored at the end */

 70:   if (format==PETSC_VIEWER_HDF5_MAT && !mat->symmetric) {
 71:     /* Swap row and columns layout for unallocated matrix. I want to avoid calling MatTranspose() just to transpose sparsity pattern and layout. */
 72:     if (!mat->preallocated) {
 73:       PetscLayout tmp;
 74:       tmp = mat->rmap; mat->rmap = mat->cmap; mat->cmap = tmp;
 75:     } else SETERRQ(comm,PETSC_ERR_SUP,"Not for preallocated matrix - we would need to transpose it here which we want to avoid");
 76:   }

 78:   /* If global sizes are set, check if they are consistent with that given in the file */
 79:   if (mat->rmap->N >= 0 && mat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Inconsistent # of rows: Matrix in file has (%D) and input matrix has (%D)",mat->rmap->N,M);
 80:   if (mat->cmap->N >= 0 && mat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Inconsistent # of cols: Matrix in file has (%D) and input matrix has (%D)",mat->cmap->N,N);

 82:   /* Determine ownership of all (block) rows and columns */
 83:   mat->rmap->N = M;
 84:   mat->cmap->N = N;
 85:   PetscLayoutSetUp(mat->rmap);
 86:   PetscLayoutSetUp(mat->cmap);
 87:   m = mat->rmap->n;

 89:   /* Read array i (array of row indices) */
 90:   PetscMalloc1(m+1, &i); /* allocate i with one more position for local number of nonzeros on each rank */
 91:   i[0] = i[m] = 0; /* make the last entry always defined - the code block below overwrites it just on last rank */
 92:   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 */
 93:   M++;
 94:   ISCreate(comm,&is_i);
 95:   PetscObjectSetName((PetscObject)is_i,i_name);
 96:   PetscLayoutSetLocalSize(is_i->map,m);
 97:   PetscLayoutSetSize(is_i->map,M);
 98:   ISLoad(is_i,viewer);
 99:   ISGetIndices(is_i,&i_glob);
100:   PetscArraycpy(i,i_glob,m);

102:   /* Reset m and M to the matrix sizes */
103:   m = mat->rmap->n;
104:   M--;

106:   /* Create PetscLayout for j and a vectors; construct ranges first */
107:   PetscMalloc1(size+1, &range);
108:   MPI_Allgather(i, 1, MPIU_INT, range, 1, MPIU_INT, comm);
109:   /* Last rank has global number of nonzeros (= length of j and a arrays) in i[m] (last i entry) so broadcast it */
110:   range[size] = i[m];
111:   MPI_Bcast(&range[size], 1, MPIU_INT, size-1, comm);
112:   for (p=size-1; p>0; p--) {
113:     if (!range[p]) range[p] = range[p+1]; /* for ranks with 0 rows, take the value from the next processor */
114:   }
115:   i[m] = range[rank+1]; /* i[m] (last i entry) is equal to next rank's offset */
116:   /* Deduce rstart, rend, n and N from the ranges */
117:   PetscLayoutCreateFromRanges(comm,range,PETSC_OWN_POINTER,1,&jmap);

119:   /* Convert global to local indexing of rows */
120:   for (p=1; p<m+1; ++p) i[p] -= i[0];
121:   i[0] = 0;

123:   /* Read array j (array of column indices) */
124:   ISCreate(comm,&is_j);
125:   PetscObjectSetName((PetscObject)is_j,j_name);
126:   PetscLayoutDuplicate(jmap,&is_j->map);
127:   ISLoad(is_j,viewer);
128:   ISGetIndices(is_j,&j);

130:   /* Read array a (array of values) */
131:   VecCreate(comm,&vec_a);
132:   PetscObjectSetName((PetscObject)vec_a,a_name);
133:   PetscLayoutDuplicate(jmap,&vec_a->map);
134:   VecLoad(vec_a,viewer);
135:   VecGetArrayRead(vec_a,&a);

137:   /* populate matrix */
138:   if (!((PetscObject)mat)->type_name) {
139:     MatSetType(mat,MATAIJ);
140:   }
141:   MatSeqAIJSetPreallocationCSR(mat,i,j,a);
142:   MatMPIAIJSetPreallocationCSR(mat,i,j,a);
143:   /*
144:   MatSeqBAIJSetPreallocationCSR(mat,bs,i,j,a);
145:   MatMPIBAIJSetPreallocationCSR(mat,bs,i,j,a);
146:   */

148:   if (format==PETSC_VIEWER_HDF5_MAT && !mat->symmetric) {
149:     /* Transpose the input matrix back */
150:     MatTranspose(mat,MAT_INPLACE_MATRIX,&mat);
151:   }

153:   PetscViewerHDF5PopGroup(viewer);
154:   PetscFree(i_name);
155:   PetscFree(j_name);
156:   PetscFree(a_name);
157:   PetscFree(c_name);
158:   PetscLayoutDestroy(&jmap);
159:   PetscFree(i);
160:   ISRestoreIndices(is_i,&i_glob);
161:   ISRestoreIndices(is_j,&j);
162:   VecRestoreArrayRead(vec_a,&a);
163:   ISDestroy(&is_i);
164:   ISDestroy(&is_j);
165:   VecDestroy(&vec_a);
166:   return(0);
167: }
168: #endif