class FullMatrix : public MatrixBase
{
public:
+
/**
* Declare type for container size.
*/
typedef types::global_dof_index size_type;
+
+ /**
+ * Default constructor. Create an empty matrix.
+ */
+ FullMatrix ();
+
+
/**
- * Create a full matrix of dimensions
- * @p m times @p n.
+ * Create a full matrix of dimensions @p m times @p n.
*/
FullMatrix (const size_type m,
const size_type n);
+
/**
- * Return a reference to the MPI
- * communicator object in use with this
- * matrix. Since this is a sequential
- * matrix, it returns the MPI_COMM_SELF
- * communicator.
+ * Throw away the present matrix and generate one that has the
+ * same properties as if it were created by the constructor of
+ * this class with the same argument list as the present function.
+ */
+ void reinit (const size_type m,
+ const size_type n);
+
+
+ /**
+ * Return a reference to the MPI communicator object in use with
+ * this matrix. Since this is a sequential matrix, it returns the
+ * MPI_COMM_SELF communicator.
*/
virtual const MPI_Comm &get_mpi_communicator () const;
+
+ private:
+
+ /**
+ * Do the actual work for the respective reinit() function and the
+ * matching constructor, i.e. create a matrix. Getting rid of the
+ * previous matrix is left to the caller.
+ */
+ void do_reinit (const size_type m,
+ const size_type n);
+
};
/*@}*/
namespace PETScWrappers
{
+
+ FullMatrix::FullMatrix ()
+ {
+ const int m=0, n=0;
+ const int ierr
+ = MatCreateSeqDense (PETSC_COMM_SELF, m, n, PETSC_NULL,
+ &matrix);
+
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+ }
+
FullMatrix::FullMatrix (const size_type m,
const size_type n)
{
+ do_reinit (m, n);
+ }
+
+ void
+ FullMatrix::reinit (const size_type m,
+ const size_type n)
+ {
+ // get rid of old matrix and generate a
+ // new one
+#if DEAL_II_PETSC_VERSION_LT(3,2,0)
+ const int ierr = MatDestroy (matrix);
+#else
+ const int ierr = MatDestroy (&matrix);
+#endif
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+ do_reinit (m, n);
+ }
+
+ void
+ FullMatrix::do_reinit (const size_type m,
+ const size_type n)
+ {
+ // use the call sequence indicating only a maximal number of
+ // elements per row for all rows globally
const int ierr
- = MatCreateSeqDense(PETSC_COMM_SELF, m, n, PETSC_NULL,
- &matrix);
+ = MatCreateSeqDense (PETSC_COMM_SELF, m, n, PETSC_NULL,
+ &matrix);
AssertThrow (ierr == 0, ExcPETScError(ierr));
}
+
const MPI_Comm &
FullMatrix::get_mpi_communicator () const
{