/**
* Destructor.
*/
- ~BlockSparseMatrix() override = default;
+ ~BlockSparseMatrix() override;
/**
* Pseudo copy operator only copying empty objects. The sizes of the
* protected.
*/
using BlockMatrixBase<SparseMatrix>::clear;
+
+ /**
+ * Conversion operator to gain access to the underlying PETSc type. If you
+ * do this, you cut this class off some information it may need, so this
+ * conversion operator should only be used if you know what you do. In
+ * particular, it should only be used for read-only operations into the
+ * matrix.
+ */
+ operator Mat() const;
+
+ /**
+ * Return a reference to the underlying PETSc type. It can be used to
+ * modify the underlying data, so use it only when you know what you
+ * are doing.
+ */
+ Mat &
+ petsc_matrix();
+
+ private:
+ Mat petsc_nest_matrix = nullptr;
};
// ---------------------------------------------------------------------
#include <deal.II/lac/petsc_block_sparse_matrix.h>
+#include <deal.II/lac/petsc_compatibility.h>
#ifdef DEAL_II_WITH_PETSC
return *this;
}
+ BlockSparseMatrix::~BlockSparseMatrix()
+ {
+ PetscErrorCode ierr = destroy_matrix(petsc_nest_matrix);
+ AssertThrow(ierr == 0, ExcPETScError(ierr));
+ }
# ifndef DOXYGEN
void
BlockSparseMatrix::collect_sizes()
{
BaseClass::collect_sizes();
+
+ auto m = this->n_block_cols();
+ auto n = this->n_block_cols();
+
+ PetscErrorCode ierr = destroy_matrix(petsc_nest_matrix);
+ AssertThrow(ierr == 0, ExcPETScError(ierr));
+ std::vector<Mat> psub_objects(m * n);
+ for (unsigned int r = 0; r < m; r++)
+ for (unsigned int c = 0; c < n; c++)
+ psub_objects[r * n + c] = this->sub_objects[r][c]->petsc_matrix();
+ ierr = MatCreateNest(get_mpi_communicator(),
+ m,
+ NULL,
+ n,
+ NULL,
+ psub_objects.data(),
+ &petsc_nest_matrix);
+ AssertThrow(ierr == 0, ExcPETScError(ierr));
}
std::vector<IndexSet>
return block(0, 0).get_mpi_communicator();
}
+ BlockSparseMatrix::operator Mat() const
+ {
+ return petsc_nest_matrix;
+ }
+
+ Mat &
+ BlockSparseMatrix::petsc_matrix()
+ {
+ return petsc_nest_matrix;
+ }
+
} // namespace MPI
} // namespace PETScWrappers
pbsm.reinit(rows, cols, bdsp, MPI_COMM_WORLD);
deallog << "nonzeros BlockSparseMatrix: " << pbsm.n_nonzero_elements()
<< std::endl;
+
+ // Extract the PETSc MATNEST and use print from PETScWrappers::MatrixBase
+ PETScWrappers::MatrixBase tmp(pbsm.petsc_matrix());
+ tmp.print(deallog.get_file_stream());
}