Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
petsc_full_matrix.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2004 - 2018 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
17
18#ifdef DEAL_II_WITH_PETSC
19
22
24
25namespace PETScWrappers
26{
28 {
29 // empty constructor generate an empty matrix
30 do_reinit(0, 0);
31 }
32
34 {
35 do_reinit(m, n);
36 }
37
38 void
40 {
41 // get rid of old matrix and generate a
42 // new one
43 const PetscErrorCode ierr = destroy_matrix(matrix);
44 AssertThrow(ierr == 0, ExcPETScError(ierr));
45
46 do_reinit(m, n);
47 }
48
49 void
51 {
52 // use the call sequence indicating only a maximal number of
53 // elements per row for all rows globally
54 const PetscErrorCode ierr =
55 MatCreateSeqDense(PETSC_COMM_SELF, m, n, nullptr, &matrix);
56 AssertThrow(ierr == 0, ExcPETScError(ierr));
57 }
58
59
60 const MPI_Comm &
62 {
63 static const MPI_Comm communicator = MPI_COMM_SELF;
64 return communicator;
65 }
66} // namespace PETScWrappers
67
68
70
71#endif // DEAL_II_WITH_PETSC
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
virtual const MPI_Comm & get_mpi_communicator() const override
void reinit(const size_type m, const size_type n)
void do_reinit(const size_type m, const size_type n)
PetscErrorCode destroy_matrix(Mat &matrix)