17 #include <deal.II/lac/block_matrix_array.h> 18 #include <deal.II/lac/block_vector.h> 19 #include <deal.II/lac/petsc_block_vector.h> 20 #include <deal.II/lac/trilinos_parallel_block_vector.h> 21 #include <deal.II/lac/vector.h> 23 DEAL_II_NAMESPACE_OPEN
26 template <
typename number,
typename BlockVectorType>
40 template <
typename number,
typename BlockVectorType>
48 template <
typename number,
typename BlockVectorType>
56 template <
typename number,
typename BlockVectorType>
58 const unsigned int n_block_rows,
59 const unsigned int n_block_cols)
60 : block_rows(n_block_rows)
61 , block_cols(n_block_cols)
65 template <
typename number,
typename BlockVectorType>
68 const unsigned int n_block_rows,
69 const unsigned int n_block_cols)
71 block_rows = n_block_rows;
72 block_cols = n_block_cols;
77 template <
typename number,
typename BlockVectorType>
80 const unsigned int n_block_rows,
81 const unsigned int n_block_cols)
84 block_rows = n_block_rows;
85 block_cols = n_block_cols;
90 template <
typename number,
typename BlockVectorType>
98 template <
typename number,
typename BlockVectorType>
101 BlockVectorType & dst,
102 const BlockVectorType &src)
const 105 Assert(dst.n_blocks() == block_rows,
107 Assert(src.n_blocks() == block_cols,
112 typename BlockVectorType::BlockType &aux = *p_aux;
114 typename std::vector<Entry>::const_iterator m = entries.begin();
115 typename std::vector<Entry>::const_iterator end = entries.end();
117 for (; m != end; ++m)
119 aux.reinit(dst.block(m->row));
121 m->matrix->Tvmult(aux, src.block(m->col));
123 m->matrix->vmult(aux, src.block(m->col));
124 dst.block(m->row).add(m->prefix, aux);
130 template <
typename number,
typename BlockVectorType>
133 BlockVectorType & dst,
134 const BlockVectorType &src)
const 142 template <
typename number,
typename BlockVectorType>
145 BlockVectorType & dst,
146 const BlockVectorType &src)
const 149 Assert(dst.n_blocks() == block_cols,
151 Assert(src.n_blocks() == block_rows,
154 typename std::vector<Entry>::const_iterator m = entries.begin();
155 typename std::vector<Entry>::const_iterator end = entries.end();
159 typename BlockVectorType::BlockType &aux = *p_aux;
161 for (; m != end; ++m)
163 aux.reinit(dst.block(m->col));
165 m->matrix->vmult(aux, src.block(m->row));
167 m->matrix->Tvmult(aux, src.block(m->row));
168 dst.block(m->col).add(m->prefix, aux);
174 template <
typename number,
typename BlockVectorType>
177 BlockVectorType & dst,
178 const BlockVectorType &src)
const 181 Tvmult_add(dst, src);
186 template <
typename number,
typename BlockVectorType>
189 const BlockVectorType &u,
190 const BlockVectorType &v)
const 193 Assert(u.n_blocks() == block_rows,
195 Assert(v.n_blocks() == block_cols,
200 typename BlockVectorType::BlockType &aux = *p_aux;
202 typename std::vector<Entry>::const_iterator m;
203 typename std::vector<Entry>::const_iterator end = entries.end();
207 for (
unsigned int i = 0; i < block_rows; ++i)
209 aux.reinit(u.block(i));
210 for (m = entries.begin(); m != end; ++m)
215 m->matrix->Tvmult_add(aux, v.block(m->col));
217 m->matrix->vmult(aux, v.block(m->col));
219 result += u.block(i) * aux;
227 template <
typename number,
typename BlockVectorType>
230 const BlockVectorType &u)
const 232 return matrix_scalar_product(u, u);
237 template <
typename number,
typename BlockVectorType>
246 template <
typename number,
typename BlockVectorType>
257 template <
typename number,
typename BlockVectorType>
264 template <
typename number,
typename BlockVectorType>
266 const unsigned int block_rows)
272 template <
typename number,
typename BlockVectorType>
280 template <
typename number,
typename BlockVectorType>
283 BlockVectorType &dst,
287 typename std::vector<
289 m = this->entries.begin();
290 typename std::vector<
292 end = this->entries.end();
293 std::vector<
typename std::vector<
299 typename BlockVectorType::BlockType &aux = *p_aux;
301 aux.reinit(dst.block(row_num),
true);
305 for (; m != end; ++m)
316 if (((j > i) && !backward) || ((j < i) && backward))
320 diagonals.push_back(m);
325 m->matrix->Tvmult(aux, dst.block(j));
327 m->matrix->vmult(aux, dst.block(j));
328 dst.block(i).add(-1.0 * m->prefix, aux);
331 Assert(diagonals.size() != 0, ExcNoDiagonal(row_num));
336 if (diagonals.size() == 1)
339 diagonals[0]->matrix->Tvmult(aux, dst.block(row_num));
341 diagonals[0]->matrix->vmult(aux, dst.block(row_num));
342 dst.block(row_num).equ(diagonals[0]->prefix, aux);
347 for (
size_type i = 0; i < diagonals.size(); ++i)
355 m->matrix->Tvmult_add(aux, dst.block(row_num));
357 m->matrix->vmult_add(aux, dst.block(row_num));
360 dst.block(row_num) = aux;
366 template <
typename number,
typename BlockVectorType>
369 BlockVectorType & dst,
370 const BlockVectorType &src)
const 372 Assert(dst.n_blocks() == n_block_rows(),
374 Assert(src.n_blocks() == n_block_cols(),
385 template <
typename number,
typename BlockVectorType>
388 BlockVectorType & dst,
389 const BlockVectorType &src)
const 391 Assert(dst.n_blocks() == n_block_rows(),
393 Assert(src.n_blocks() == n_block_cols(),
400 for (
unsigned int i = n_block_rows(); i > 0;)
405 for (
unsigned int i = 0; i < n_block_rows(); ++i)
410 template <
typename number,
typename BlockVectorType>
414 const BlockVectorType &)
const 420 template <
typename number,
typename BlockVectorType>
424 const BlockVectorType &)
const 434 #ifdef DEAL_II_WITH_TRILINOS 443 #ifdef DEAL_II_WITH_PETSC 444 # ifdef PETSC_USE_COMPLEX 463 DEAL_II_NAMESPACE_CLOSE
number matrix_norm_square(const BlockVectorType &u) const
PointerMatrixBase< typename BlockVectorType::BlockType > * matrix
void Tvmult(BlockVectorType &dst, const BlockVectorType &src) const
void vmult(BlockVectorType &dst, const BlockVectorType &src) const
types::global_dof_index size_type
void vmult_add(BlockVectorType &dst, const BlockVectorType &src) const
unsigned int n_block_cols() const
void initialize(const unsigned int n_block_rows, const unsigned int n_block_cols)
#define Assert(cond, exc)
void vmult_add(BlockVectorType &dst, const BlockVectorType &src) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int n_block_rows() const
void Tvmult(BlockVectorType &dst, const BlockVectorType &src) const
Entry(const MatrixType &matrix, size_type row, size_type col, number prefix, bool transpose)
void reinit(const unsigned int n_block_rows, const unsigned int n_block_cols)
BlockTrianglePrecondition()
void Tvmult_add(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_add(BlockVectorType &dst, const BlockVectorType &src) const
static ::ExceptionBase & ExcNotImplemented()
void reinit(const unsigned int n_block_rows)
void vmult(BlockVectorType &dst, const BlockVectorType &src) const
void do_row(BlockVectorType &dst, size_type row_num) const
number matrix_scalar_product(const BlockVectorType &u, const BlockVectorType &v) const