xref: /petsc/src/mat/utils/matio.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
13ea6fe3dSLisandro Dalcin #include <petscviewer.h>
23ea6fe3dSLisandro Dalcin #include <petsc/private/matimpl.h>
33ea6fe3dSLisandro Dalcin 
4*9371c9d4SSatish Balay PetscErrorCode MatView_Binary_BlockSizes(Mat mat, PetscViewer viewer) {
53ea6fe3dSLisandro Dalcin   FILE       *info;
63ea6fe3dSLisandro Dalcin   PetscMPIInt rank;
73ea6fe3dSLisandro Dalcin   PetscInt    rbs, cbs;
83ea6fe3dSLisandro Dalcin 
93ea6fe3dSLisandro Dalcin   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(mat, &rbs, &cbs));
119566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetInfoPointer(viewer, &info));
129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
13dd400576SPatrick Sanan   if (rank == 0 && info) {
149566063dSJacob Faibussowitsch     if (rbs != cbs) PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "-matload_block_size %" PetscInt_FMT ",%" PetscInt_FMT "\n", rbs, cbs));
159566063dSJacob Faibussowitsch     else PetscCall(PetscFPrintf(PETSC_COMM_SELF, info, "-matload_block_size %" PetscInt_FMT "\n", rbs));
163ea6fe3dSLisandro Dalcin   }
173ea6fe3dSLisandro Dalcin   PetscFunctionReturn(0);
183ea6fe3dSLisandro Dalcin }
193ea6fe3dSLisandro Dalcin 
20*9371c9d4SSatish Balay PetscErrorCode MatLoad_Binary_BlockSizes(Mat mat, PetscViewer viewer) {
213ea6fe3dSLisandro Dalcin   PetscInt  rbs, cbs, bs[2], n = 2;
223ea6fe3dSLisandro Dalcin   PetscBool set;
233ea6fe3dSLisandro Dalcin 
243ea6fe3dSLisandro Dalcin   PetscFunctionBegin;
253ea6fe3dSLisandro Dalcin   /* get current block sizes */
269566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(mat, &rbs, &cbs));
27*9371c9d4SSatish Balay   bs[0] = rbs;
28*9371c9d4SSatish Balay   bs[1] = cbs;
293ea6fe3dSLisandro Dalcin   /* get block sizes from the options database */
30d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)viewer), NULL, "Options for loading matrix block size", "Mat");
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-matload_block_size", "Set the block size used to store the matrix", "MatLoad", bs, &n, &set));
32d0609cedSBarry Smith   PetscOptionsEnd();
333ea6fe3dSLisandro Dalcin   if (!set) PetscFunctionReturn(0);
343ea6fe3dSLisandro Dalcin   if (n == 1) bs[1] = bs[0]; /* to support -matload_block_size <bs> */
353ea6fe3dSLisandro Dalcin   /* set matrix block sizes */
363ea6fe3dSLisandro Dalcin   if (bs[0] > 0) rbs = bs[0];
373ea6fe3dSLisandro Dalcin   if (bs[1] > 0) cbs = bs[1];
389566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(mat, rbs, cbs));
393ea6fe3dSLisandro Dalcin   PetscFunctionReturn(0);
403ea6fe3dSLisandro Dalcin }
41