Reference documentation for deal.II version 8.5.1
trilinos_block_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2017 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__trilinos_block_sparse_matrix_h
17 #define dealii__trilinos_block_sparse_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/base/table.h>
25 # include <deal.II/base/template_constraints.h>
26 # include <deal.II/lac/block_matrix_base.h>
27 # include <deal.II/lac/trilinos_sparse_matrix.h>
28 # include <deal.II/lac/trilinos_block_vector.h>
29 # include <deal.II/lac/full_matrix.h>
30 # include <deal.II/lac/exceptions.h>
31 
32 # include <cmath>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 // forward declarations
38 template <typename number> class BlockSparseMatrix;
39 
40 
41 namespace TrilinosWrappers
42 {
43 
71  class BlockSparseMatrix : public BlockMatrixBase<SparseMatrix>
72  {
73  public:
78 
83 
88  typedef BaseClass::pointer pointer;
89  typedef BaseClass::const_pointer const_pointer;
90  typedef BaseClass::reference reference;
91  typedef BaseClass::const_reference const_reference;
92  typedef BaseClass::size_type size_type;
95 
108 
113 
120 
131  operator = (const double d);
132 
146  void reinit (const size_type n_block_rows,
147  const size_type n_block_columns);
148 
154  template <typename BlockSparsityPatternType>
155  void reinit (const std::vector<Epetra_Map> &input_maps,
156  const BlockSparsityPatternType &block_sparsity_pattern,
157  const bool exchange_data = false);
158 
164  template <typename BlockSparsityPatternType>
165  void reinit (const std::vector<IndexSet> &input_maps,
166  const BlockSparsityPatternType &block_sparsity_pattern,
167  const MPI_Comm &communicator = MPI_COMM_WORLD,
168  const bool exchange_data = false);
169 
175  template <typename BlockSparsityPatternType>
176  void reinit (const BlockSparsityPatternType &block_sparsity_pattern);
177 
186  void reinit (const std::vector<Epetra_Map> &input_maps,
187  const ::BlockSparseMatrix<double> &deal_ii_sparse_matrix,
188  const double drop_tolerance=1e-13) DEAL_II_DEPRECATED;
189 
197  void reinit (const ::BlockSparseMatrix<double> &deal_ii_sparse_matrix,
198  const double drop_tolerance=1e-13);
199 
207  bool is_compressed () const;
208 
219  void collect_sizes ();
220 
225  size_type n_nonzero_elements () const;
226 
236  std::vector<Epetra_Map> domain_partitioner () const DEAL_II_DEPRECATED;
237 
247  std::vector<Epetra_Map> range_partitioner () const DEAL_II_DEPRECATED;
248 
249 
256  template <typename VectorType1, typename VectorType2>
257  void vmult (VectorType1 &dst,
258  const VectorType2 &src) const;
259 
265  template <typename VectorType1, typename VectorType2>
266  void Tvmult (VectorType1 &dst,
267  const VectorType2 &src) const;
268 
281  TrilinosScalar residual (MPI::BlockVector &dst,
282  const MPI::BlockVector &x,
283  const MPI::BlockVector &b) const;
284 
299  TrilinosScalar residual (BlockVector &dst,
300  const BlockVector &x,
301  const BlockVector &b) const;
302 
310  TrilinosScalar residual (MPI::BlockVector &dst,
311  const MPI::Vector &x,
312  const MPI::BlockVector &b) const;
313 
321  TrilinosScalar residual (BlockVector &dst,
322  const Vector &x,
323  const BlockVector &b) const;
324 
332  TrilinosScalar residual (MPI::Vector &dst,
333  const MPI::BlockVector &x,
334  const MPI::Vector &b) const;
335 
343  TrilinosScalar residual (Vector &dst,
344  const BlockVector &x,
345  const Vector &b) const;
346 
354  TrilinosScalar residual (VectorBase &dst,
355  const VectorBase &x,
356  const VectorBase &b) const;
357 
363 
373  int, int, int, int,
374  << "The blocks [" << arg1 << ',' << arg2 << "] and ["
375  << arg3 << ',' << arg4 << "] have differing row numbers.");
376 
381  int, int, int, int,
382  << "The blocks [" << arg1 << ',' << arg2 << "] and ["
383  << arg3 << ',' << arg4 << "] have differing column numbers.");
385 
386  private:
390  template <typename VectorType1, typename VectorType2>
391  void vmult (VectorType1 &dst,
392  const VectorType2 &src,
393  const bool transpose,
394  const ::internal::bool2type<true>,
395  const ::internal::bool2type<true>) const;
396 
401  template <typename VectorType1, typename VectorType2>
402  void vmult (VectorType1 &dst,
403  const VectorType2 &src,
404  const bool transpose,
405  const ::internal::bool2type<false>,
406  const ::internal::bool2type<true>) const;
407 
412  template <typename VectorType1, typename VectorType2>
413  void vmult (VectorType1 &dst,
414  const VectorType2 &src,
415  const bool transpose,
416  const ::internal::bool2type<true>,
417  const ::internal::bool2type<false>) const;
418 
424  template <typename VectorType1, typename VectorType2>
425  void vmult (VectorType1 &dst,
426  const VectorType2 &src,
427  const bool transpose,
428  const ::internal::bool2type<false>,
429  const ::internal::bool2type<false>) const;
430  };
431 
432 
433 
436 // ------------- inline and template functions -----------------
437 
438 
439 
440  inline
442  BlockSparseMatrix::operator = (const double d)
443  {
445 
446  for (size_type r=0; r<this->n_block_rows(); ++r)
447  for (size_type c=0; c<this->n_block_cols(); ++c)
448  this->block(r,c) = d;
449 
450  return *this;
451  }
452 
453 
454 
455  inline
456  bool
458  {
459  bool compressed = true;
460  for (size_type row=0; row<n_block_rows(); ++row)
461  for (size_type col=0; col<n_block_cols(); ++col)
462  if (block(row, col).is_compressed() == false)
463  {
464  compressed = false;
465  break;
466  }
467 
468  return compressed;
469  }
470 
471 
472 
473  template <typename VectorType1, typename VectorType2>
474  inline
475  void
476  BlockSparseMatrix::vmult (VectorType1 &dst,
477  const VectorType2 &src) const
478  {
479  vmult(dst, src, false,
482  }
483 
484 
485 
486  template <typename VectorType1, typename VectorType2>
487  inline
488  void
489  BlockSparseMatrix::Tvmult (VectorType1 &dst,
490  const VectorType2 &src) const
491  {
492  vmult(dst, src, true,
495  }
496 
497 
498 
499  template <typename VectorType1, typename VectorType2>
500  inline
501  void
502  BlockSparseMatrix::vmult (VectorType1 &dst,
503  const VectorType2 &src,
504  const bool transpose,
507  {
508  if (transpose == true)
510  else
511  BaseClass::vmult_block_block (dst, src);
512  }
513 
514 
515 
516 
517  template <typename VectorType1, typename VectorType2>
518  inline
519  void
520  BlockSparseMatrix::vmult (VectorType1 &dst,
521  const VectorType2 &src,
522  const bool transpose,
525  {
526  if (transpose == true)
528  else
530  }
531 
532 
533 
534  template <typename VectorType1, typename VectorType2>
535  inline
536  void
537  BlockSparseMatrix::vmult (VectorType1 &dst,
538  const VectorType2 &src,
539  const bool transpose,
542  {
543  if (transpose == true)
545  else
547  }
548 
549 
550 
551  template <typename VectorType1, typename VectorType2>
552  inline
553  void
554  BlockSparseMatrix::vmult (VectorType1 &dst,
555  const VectorType2 &src,
556  const bool transpose,
559  {
560  if (transpose == true)
562  else
564  }
565 
566 
567 #ifdef DEAL_II_WITH_CXX11
568 
569 
570  namespace internal
571  {
572  namespace BlockLinearOperator
573  {
574 
590  template<typename PayloadBlockType>
592  {
593  public:
597  typedef PayloadBlockType BlockType;
598 
608  template <typename... Args>
609  TrilinosBlockPayload (const Args &...)
610  {
611  static_assert(typeid(PayloadBlockType)==typeid(internal::LinearOperator::TrilinosPayload),
612  "TrilinosBlockPayload can only accept a payload of type TrilinosPayload.");
613  }
614  };
615 
616  } /*namespace BlockLinearOperator*/
617  } /* namespace internal */
618 
619 #endif // DEAL_II_WITH_CXX11
620 
621 }/* namespace TrilinosWrappers */
622 
623 
624 DEAL_II_NAMESPACE_CLOSE
625 
626 #endif // DEAL_II_WITH_TRILINOS
627 
628 #endif // dealii__trilinos_block_sparse_matrix_h
void vmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcIncompatibleRowNumbers(int arg1, int arg2, int arg3, int arg4)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
BlockMatrixBase< SparseMatrix > BaseClass
void vmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
STL namespace.
void Tvmult_block_nonblock(BlockVectorType &dst, const VectorType &src) const
void Tvmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_block(VectorType &dst, const BlockVectorType &src) const
void vmult_block_block(BlockVectorType &dst, const BlockVectorType &src) const
void Tvmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcIncompatibleColNumbers(int arg1, int arg2, int arg3, int arg4)
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
unsigned int n_block_cols() const
void vmult(VectorType1 &dst, const VectorType2 &src) const
#define Assert(cond, exc)
Definition: exceptions.h:313
std::vector< Epetra_Map > range_partitioner() const 1
void reinit(const size_type n_block_rows, const size_type n_block_columns)
BlockType::value_type value_type
void Tvmult(VectorType1 &dst, const VectorType2 &src) const
void vmult_nonblock_nonblock(VectorType &dst, const VectorType &src) const
BlockType & block(const unsigned int row, const unsigned int column)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
std::vector< Epetra_Map > domain_partitioner() const 1
unsigned int n_block_rows() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const