17 #include <deal.II/lac/block_matrix_array.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/lac/block_vector.h> 20 #include <deal.II/lac/trilinos_parallel_block_vector.h> 21 #include <deal.II/lac/trilinos_parallel_block_vector.h> 23 #include <deal.II/lac/petsc_parallel_block_vector.h> 25 DEAL_II_NAMESPACE_OPEN
28 template <
typename number,
typename BlockVectorType>
43 template <
typename number,
typename BlockVectorType>
51 template <
typename number,
typename BlockVectorType>
59 template <
typename number,
typename BlockVectorType>
61 (
const unsigned int n_block_rows,
62 const unsigned int n_block_cols)
63 : block_rows (n_block_rows),
64 block_cols (n_block_cols)
68 template <
typename number,
typename BlockVectorType>
71 (
const unsigned int n_block_rows,
72 const unsigned int n_block_cols)
74 block_rows = n_block_rows;
75 block_cols = n_block_cols;
80 template <
typename number,
typename BlockVectorType>
83 (
const unsigned int n_block_rows,
84 const unsigned int n_block_cols)
87 block_rows = n_block_rows;
88 block_cols = n_block_cols;
93 template <
typename number,
typename BlockVectorType>
101 template <
typename number,
typename BlockVectorType>
104 const BlockVectorType &src)
const 107 Assert (dst.n_blocks() == block_rows,
109 Assert (src.n_blocks() == block_cols,
113 typename BlockVectorType::BlockType &aux = *p_aux;
115 typename std::vector<Entry>::const_iterator m = entries.begin();
116 typename std::vector<Entry>::const_iterator end = entries.end();
118 for (; m != end ; ++m)
120 aux.reinit(dst.block(m->row));
122 m->matrix->Tvmult(aux, src.block(m->col));
124 m->matrix->vmult(aux, src.block(m->col));
125 dst.block(m->row).add (m->prefix, aux);
132 template <
typename number,
typename BlockVectorType>
135 const BlockVectorType &src)
const 138 vmult_add (dst, src);
144 template <
typename number,
typename BlockVectorType>
147 const BlockVectorType &src)
const 150 Assert (dst.n_blocks() == block_cols,
152 Assert (src.n_blocks() == block_rows,
155 typename std::vector<Entry>::const_iterator m = entries.begin();
156 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 const BlockVectorType &src)
const 180 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,
199 typename BlockVectorType::BlockType &aux = *p_aux;
201 typename std::vector<Entry>::const_iterator m;
202 typename std::vector<Entry>::const_iterator end = entries.end();
206 for (
unsigned int i=0; i<block_rows; ++i)
208 aux.reinit(u.block(i));
209 for (m = entries.begin(); m != end ; ++m)
214 m->matrix->Tvmult_add(aux, v.block(m->col));
216 m->matrix->vmult(aux, v.block(m->col));
218 result += u.block(i)*aux;
226 template <
typename number,
typename BlockVectorType>
229 (
const BlockVectorType &u)
const 231 return matrix_scalar_product(u,u);
236 template <
typename number,
typename BlockVectorType>
245 template <
typename number,
typename BlockVectorType>
256 template <
typename number,
typename BlockVectorType>
263 template <
typename number,
typename BlockVectorType>
265 (
const unsigned int block_rows)
272 template <
typename number,
typename BlockVectorType>
275 (
const unsigned int n)
281 template <
typename number,
typename BlockVectorType>
284 (BlockVectorType &dst,
288 typename std::vector<typename BlockMatrixArray<number,BlockVectorType>::Entry>::const_iterator
289 m = this->entries.begin();
290 typename std::vector<typename BlockMatrixArray<number,BlockVectorType>::Entry>::const_iterator
291 end = this->entries.end();
292 std::vector<typename std::vector<typename BlockMatrixArray<number,BlockVectorType>::Entry>::const_iterator>
296 typename BlockVectorType::BlockType &aux = *p_aux;
298 aux.reinit(dst.block(row_num),
true);
302 for (; m != end ; ++m)
313 if (((j > i) && !backward) || ((j < i) && backward))
317 diagonals.push_back(m);
322 m->matrix->Tvmult(aux, dst.block(j));
324 m->matrix->vmult(aux, dst.block(j));
325 dst.block(i).add (-1.0 * m->prefix, aux);
328 Assert (diagonals.size() != 0, ExcNoDiagonal(row_num));
333 if (diagonals.size() == 1)
336 diagonals[0]->matrix->Tvmult(aux, dst.block(row_num));
338 diagonals[0]->matrix->vmult(aux, dst.block(row_num));
339 dst.block(row_num).equ (diagonals[0]->prefix, aux);
344 for (
size_type i=0; i<diagonals.size(); ++i)
352 m->matrix->Tvmult_add(aux, dst.block(row_num));
354 m->matrix->vmult_add(aux, dst.block(row_num));
357 dst.block(row_num) = aux;
363 template <
typename number,
typename BlockVectorType>
366 (BlockVectorType &dst,
367 const BlockVectorType &src)
const 369 Assert (dst.n_blocks() == n_block_rows(),
371 Assert (src.n_blocks() == n_block_cols(),
382 template <
typename number,
typename BlockVectorType>
385 (BlockVectorType &dst,
386 const BlockVectorType &src)
const 388 Assert (dst.n_blocks() == n_block_rows(),
390 Assert (src.n_blocks() == n_block_cols(),
397 for (
unsigned int i=n_block_rows(); i>0;)
402 for (
unsigned int i=0; i<n_block_rows(); ++i)
408 template <
typename number,
typename BlockVectorType>
412 const BlockVectorType &)
const 418 template <
typename number,
typename BlockVectorType>
422 const BlockVectorType &)
const 432 #ifdef DEAL_II_WITH_TRILINOS 439 #ifdef DEAL_II_WITH_PETSC 440 #ifdef PETSC_USE_COMPLEX 453 DEAL_II_NAMESPACE_CLOSE
number matrix_norm_square(const BlockVectorType &u) const
types::global_dof_index size_type
void Tvmult(BlockVectorType &dst, const BlockVectorType &src) const
void vmult(BlockVectorType &dst, const BlockVectorType &src) const
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
PointerMatrixBase< typename BlockVectorType::BlockType > * matrix
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