16 #ifndef dealii__block_linear_operator_h 17 #define dealii__block_linear_operator_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/exceptions.h> 22 #include <deal.II/lac/linear_operator.h> 24 #ifdef DEAL_II_WITH_CXX11 26 DEAL_II_NAMESPACE_OPEN
32 template <
typename Range = BlockVector<
double>,
33 typename Domain = Range>
36 template <
typename Range = BlockVector<
double>,
37 typename Domain = Range,
38 typename BlockMatrixType>
42 template <
size_t m,
size_t n,
44 typename Domain = Range>
50 typename Domain = Range>
56 typename Domain = Range>
67 template <
typename Range = BlockVector<
double>,
68 typename Domain = Range,
69 typename BlockMatrixType>
73 template <
typename Range = BlockVector<
double>,
74 typename Domain = Range>
79 template <
typename Range = BlockVector<
double>,
80 typename Domain = Range>
128 template <
typename Range,
typename Domain>
148 Assert(
false,
ExcMessage(
"Uninitialized BlockLinearOperator<Range, Domain>::n_block_rows called"));
154 Assert(
false,
ExcMessage(
"Uninitialized BlockLinearOperator<Range, Domain>::n_block_cols called"));
158 block = [](
unsigned int,
unsigned int) -> BlockType
160 Assert(
false,
ExcMessage(
"Uninitialized BlockLinearOperator<Range, Domain>::block called"));
176 template<
typename Op>
179 *
this = block_operator<Range, Domain, Op>(op);
187 template<
size_t m,
size_t n>
190 *
this = block_operator<m, n, Range, Domain>(ops);
201 *
this = block_diagonal_operator<m, Range, Domain>(ops);
214 template <
typename Op>
217 *
this = block_operator<Range, Domain, Op>(op);
226 template <
size_t m,
size_t n>
228 operator=(
const std::array<std::array<BlockType, n>, m> &ops)
230 *
this = block_operator<m, n, Range, Domain>(ops);
243 *
this = block_diagonal_operator<m, Range, Domain>(ops);
264 std::function<BlockType(unsigned int, unsigned int)>
block;
275 template <
typename Range,
typename Domain>
277 populate_linear_operator_functions(
288 for (
unsigned int i = 0; i < m; ++i)
289 op.
block(i, 0).reinit_range_vector(v.block(i), omit_zeroing_entries);
302 for (
unsigned int i = 0; i < n; ++i)
303 op.
block(0, i).reinit_domain_vector(v.block(i), omit_zeroing_entries);
308 op.
vmult = [=](Range &v,
const Domain &u)
312 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
313 Assert(u.n_blocks() == n, ExcDimensionMismatch(u.n_blocks(), n));
315 for (
unsigned int i = 0; i < m; ++i)
317 op.
block(i, 0).vmult(v.block(i), u.block(0));
318 for (
unsigned int j = 1; j < n; ++j)
319 op.
block(i, j).vmult_add(v.block(i), u.block(j));
323 op.
vmult_add = [=](Range &v,
const Domain &u)
327 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
328 Assert(u.n_blocks() == n, ExcDimensionMismatch(u.n_blocks(), n));
330 for (
unsigned int i = 0; i < m; ++i)
331 for (
unsigned int j = 0; j < n; ++j)
332 op.
block(i, j).vmult_add(v.block(i), u.block(j));
335 op.
Tvmult = [=](Domain &v,
const Range &u)
339 Assert(v.n_blocks() == n, ExcDimensionMismatch(v.n_blocks(), n));
340 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
342 for (
unsigned int i = 0; i < n; ++i)
344 op.
block(0, i).Tvmult(v.block(i), u.block(0));
345 for (
unsigned int j = 1; j < m; ++j)
346 op.
block(j, i).Tvmult_add(v.block(i), u.block(j));
350 op.
Tvmult_add = [=](Domain &v,
const Range &u)
354 Assert(v.n_blocks() == n, ExcDimensionMismatch(v.n_blocks(), n));
355 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
357 for (
unsigned int i = 0; i < n; ++i)
358 for (
unsigned int j = 0; j < m; ++j)
359 op.
block(j, i).Tvmult_add(v.block(i), u.block(j));
383 template <
typename Range,
385 typename BlockMatrixType>
393 return_op.
n_block_rows = [&block_matrix]() ->
unsigned int 395 return block_matrix.n_block_rows();
398 return_op.
n_block_cols = [&block_matrix]() ->
unsigned int 400 return block_matrix.n_block_cols();
403 return_op.
block = [&block_matrix](
unsigned int i,
unsigned int j) -> BlockType
406 const unsigned int m = block_matrix.n_block_rows();
407 const unsigned int n = block_matrix.n_block_cols();
408 Assert(i < m, ExcIndexRange (i, 0, m));
409 Assert(j < n, ExcIndexRange (j, 0, n));
412 return BlockType(block_matrix.block(i, j));
415 internal::BlockLinearOperator::populate_linear_operator_functions(return_op);
448 template <
size_t m,
size_t n,
typename Range,
typename Domain>
452 static_assert(m > 0 && n > 0,
453 "a blocked LinearOperator must consist of at least one block");
469 return_op.
block = [ops](
unsigned int i,
unsigned int j) -> BlockType
471 Assert(i < m, ExcIndexRange (i, 0, m));
472 Assert(j < n, ExcIndexRange (j, 0, n));
477 internal::BlockLinearOperator::populate_linear_operator_functions(return_op);
498 template <
typename Range,
500 typename BlockMatrixType>
508 return_op.
n_block_rows = [&block_matrix]() ->
unsigned int 510 return block_matrix.n_block_rows();
513 return_op.
n_block_cols = [&block_matrix]() ->
unsigned int 515 return block_matrix.n_block_cols();
518 return_op.
block = [&block_matrix](
unsigned int i,
unsigned int j) -> BlockType
521 const unsigned int m = block_matrix.n_block_rows();
522 const unsigned int n = block_matrix.n_block_cols();
523 Assert(m == n, ExcDimensionMismatch(m, n));
524 Assert(i < m, ExcIndexRange (i, 0, m));
525 Assert(j < n, ExcIndexRange (j, 0, n));
528 return BlockType(block_matrix.block(i, j));
533 internal::BlockLinearOperator::populate_linear_operator_functions(return_op);
556 template <
size_t m,
typename Range,
typename Domain>
561 "a blockdiagonal LinearOperator must consist of at least one block");
565 std::array<std::array<BlockType, m>, m> new_ops;
571 for (
unsigned int i = 0; i < m; ++i)
572 for (
unsigned int j = 0; j < m; ++j)
576 new_ops[i][j] = ops[i];
583 new_ops[i][j].reinit_domain_vector = ops[j].reinit_domain_vector;
586 return block_operator<m,m,Range,Domain>(new_ops);
600 template <
size_t m,
typename Range,
typename Domain>
605 "a blockdiagonal LinearOperator must consist of at least " 609 std::array<BlockType, m> new_ops;
658 template <
typename Range,
typename Domain>
672 ExcDimensionMismatch(block_operator.
n_block_cols(), m));
674 ExcDimensionMismatch(diagonal_inverse.
n_block_rows(), m));
676 ExcDimensionMismatch(diagonal_inverse.
n_block_cols(), m));
677 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
678 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
683 diagonal_inverse.
block(0, 0).vmult(v.block(0), u.block(0));
684 for (
unsigned int i = 1; i < m; ++i)
686 auto &dst = v.block(i);
689 for (
unsigned int j = 0; j < i; ++j)
690 block_operator.
block(i, j).vmult_add(dst, v.block(j));
692 diagonal_inverse.
block(i, i).vmult(dst, dst);
700 ExcDimensionMismatch(block_operator.
n_block_cols(), m));
702 ExcDimensionMismatch(diagonal_inverse.
n_block_rows(), m));
704 ExcDimensionMismatch(diagonal_inverse.
n_block_cols(), m));
705 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
706 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
712 typename Range::BlockType *tmp = vector_memory.
alloc();
714 diagonal_inverse.
block(0, 0).vmult_add(v.block(0), u.block(0));
716 for (
unsigned int i = 1; i < m; ++i)
718 diagonal_inverse.
block(i, i).reinit_range_vector(*tmp,
true);
721 for (
unsigned int j = 0; j < i; ++j)
722 block_operator.
block(i, j).vmult_add(*tmp, v.block(j));
724 diagonal_inverse.
block(i, i).vmult_add(v.block(i),*tmp);
727 vector_memory.
free(tmp);
770 template <
typename Range,
typename Domain>
784 ExcDimensionMismatch(block_operator.
n_block_cols(), m));
786 ExcDimensionMismatch(diagonal_inverse.
n_block_rows(), m));
788 ExcDimensionMismatch(diagonal_inverse.
n_block_cols(), m));
789 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
790 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
795 diagonal_inverse.
block(m-1, m-1).vmult(v.block(m-1),u.block(m-1));
797 for (
int i = m - 2; i >= 0; --i)
799 auto &dst = v.block(i);
802 for (
unsigned int j = i + 1; j < m; ++j)
803 block_operator.
block(i, j).vmult_add(dst, v.block(j));
805 diagonal_inverse.
block(i, i).vmult(dst, dst);
813 ExcDimensionMismatch(block_operator.
n_block_cols(), m));
815 ExcDimensionMismatch(diagonal_inverse.
n_block_rows(), m));
817 ExcDimensionMismatch(diagonal_inverse.
n_block_cols(), m));
818 Assert(v.n_blocks() == m, ExcDimensionMismatch(v.n_blocks(), m));
819 Assert(u.n_blocks() == m, ExcDimensionMismatch(u.n_blocks(), m));
821 typename Range::BlockType *tmp = vector_memory.
alloc();
826 diagonal_inverse.
block(m-1, m-1).vmult_add(v.block(m-1),u.block(m-1));
828 for (
int i = m - 2; i >= 0; --i)
830 diagonal_inverse.
block(i, i).reinit_range_vector(*tmp,
true);
833 for (
unsigned int j = i + 1; j < m; ++j)
834 block_operator.
block(i, j).vmult_add(*tmp,v.block(j));
836 diagonal_inverse.
block(i, i).vmult_add(v.block(i),*tmp);
839 vector_memory.
free(tmp);
847 DEAL_II_NAMESPACE_CLOSE
849 #endif // DEAL_II_WITH_CXX11 std::function< unsigned int()> n_block_cols
BlockLinearOperator< Range, Domain > block_diagonal_operator(const BlockMatrixType &block_matrix)
::ExceptionBase & ExcMessage(std::string arg1)
BlockLinearOperator< Range, Domain > & operator=(const std::array< BlockType, m > &ops)
BlockLinearOperator< Range, Domain > block_diagonal_operator(const LinearOperator< typename Range::BlockType, typename Domain::BlockType > &op)
std::function< unsigned int()> n_block_rows
LinearOperator< Domain, Range > block_forward_substitution(const BlockLinearOperator< Range, Domain > &block_operator, const BlockLinearOperator< Domain, Range > &diagonal_inverse)
std::function< void(Range &v, bool omit_zeroing_entries)> reinit_range_vector
std::function< BlockType(unsigned int, unsigned int)> block
BlockLinearOperator(const Op &op)
BlockLinearOperator< Range, Domain > block_diagonal_operator(const std::array< LinearOperator< typename Range::BlockType, typename Domain::BlockType >, m > &ops)
std::function< void(Domain &v, const Range &u)> Tvmult_add
LinearOperator< Range, Domain > null_operator(const LinearOperator< Range, Domain > &op)
BlockLinearOperator(const std::array< BlockType, m > &ops)
virtual VectorType * alloc()
#define Assert(cond, exc)
virtual void free(const VectorType *const)
std::function< void(Domain &v, bool omit_zeroing_entries)> reinit_domain_vector
BlockLinearOperator< Range, Domain > block_operator(const std::array< std::array< LinearOperator< typename Range::BlockType, typename Domain::BlockType >, n >, m > &ops)
BlockLinearOperator< Range, Domain > block_operator(const BlockMatrixType &block_matrix)
BlockLinearOperator< Range, Domain > & operator=(const Op &op)
BlockLinearOperator< Range, Domain > & operator=(const std::array< std::array< BlockType, n >, m > &ops)
std::function< void(Range &v, const Domain &u)> vmult
LinearOperator< Domain, Range > block_back_substitution(const BlockLinearOperator< Range, Domain > &block_operator, const BlockLinearOperator< Domain, Range > &diagonal_inverse)
std::function< void(Range &v, const Domain &u)> vmult_add
BlockLinearOperator< Range, Domain > & operator=(const BlockLinearOperator< Range, Domain > &)=default
std::function< void(Domain &v, const Range &u)> Tvmult
BlockLinearOperator(const std::array< std::array< BlockType, n >, m > &ops)