Reference documentation for deal.II version 8.5.1
trilinos_block_sparse_matrix.cc
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 #include <deal.II/lac/trilinos_block_sparse_matrix.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/lac/block_sparse_matrix.h>
21 # include <deal.II/lac/block_sparsity_pattern.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 namespace TrilinosWrappers
26 {
28  {}
29 
30 
31 
33  {
34  // delete previous content of
35  // the subobjects array
36  try
37  {
38  clear ();
39  }
40  catch (...)
41  {}
42  }
43 
44 
45 
48  {
50 
51  return *this;
52  }
53 
54 
55 
56  void
58  reinit (const size_type n_block_rows,
59  const size_type n_block_columns)
60  {
61  // first delete previous content of
62  // the subobjects array
63  clear ();
64 
65  // then resize. set sizes of blocks to
66  // zero. user will later have to call
67  // collect_sizes for this
68  this->sub_objects.reinit (n_block_rows,
69  n_block_columns);
70  this->row_block_indices.reinit (n_block_rows, 0);
71  this->column_block_indices.reinit (n_block_columns, 0);
72 
73  // and reinitialize the blocks
74  for (size_type r=0; r<this->n_block_rows(); ++r)
75  for (size_type c=0; c<this->n_block_cols(); ++c)
76  {
77  BlockType *p = new BlockType();
78 
79  Assert (this->sub_objects[r][c] == 0,
81  this->sub_objects[r][c] = p;
82  }
83  }
84 
85 
86 
87  template <typename BlockSparsityPatternType>
88  void
90  reinit (const std::vector<Epetra_Map> &parallel_partitioning,
91  const BlockSparsityPatternType &block_sparsity_pattern,
92  const bool exchange_data)
93  {
94  Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_rows(),
95  ExcDimensionMismatch (parallel_partitioning.size(),
96  block_sparsity_pattern.n_block_rows()));
97  Assert (parallel_partitioning.size() == block_sparsity_pattern.n_block_cols(),
98  ExcDimensionMismatch (parallel_partitioning.size(),
99  block_sparsity_pattern.n_block_cols()));
100 
101  const size_type n_block_rows = parallel_partitioning.size();
102  (void)n_block_rows;
103 
104  Assert (n_block_rows == block_sparsity_pattern.n_block_rows(),
106  block_sparsity_pattern.n_block_rows()));
107  Assert (n_block_rows == block_sparsity_pattern.n_block_cols(),
109  block_sparsity_pattern.n_block_cols()));
110 
111 
112  // Call the other basic reinit function, ...
113  reinit (block_sparsity_pattern.n_block_rows(),
114  block_sparsity_pattern.n_block_cols());
115 
116  // ... set the correct sizes, ...
117  this->row_block_indices = block_sparsity_pattern.get_row_indices();
118  this->column_block_indices = block_sparsity_pattern.get_column_indices();
119 
120  // ... and then assign the correct
121  // data to the blocks.
122  for (size_type r=0; r<this->n_block_rows(); ++r)
123  for (size_type c=0; c<this->n_block_cols(); ++c)
124  {
125  this->sub_objects[r][c]->reinit (parallel_partitioning[r],
126  parallel_partitioning[c],
127  block_sparsity_pattern.block(r,c),
128  exchange_data);
129  }
130  }
131 
132 
133 
134  template <typename BlockSparsityPatternType>
135  void
137  reinit (const std::vector<IndexSet> &parallel_partitioning,
138  const BlockSparsityPatternType &block_sparsity_pattern,
139  const MPI_Comm &communicator,
140  const bool exchange_data)
141  {
142  std::vector<Epetra_Map> epetra_maps;
143  for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
144  epetra_maps.push_back
145  (parallel_partitioning[i].make_trilinos_map(communicator, false));
146 
147  reinit (epetra_maps, block_sparsity_pattern, exchange_data);
148 
149  }
150 
151 
152 
153  template <typename BlockSparsityPatternType>
154  void
156  reinit (const BlockSparsityPatternType &block_sparsity_pattern)
157  {
158  std::vector<Epetra_Map> parallel_partitioning;
159  for (size_type i=0; i<block_sparsity_pattern.n_block_rows(); ++i)
160  parallel_partitioning.push_back
161  (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(block_sparsity_pattern.block(i,0).n_rows()),
162  0,
164 
165  reinit (parallel_partitioning, block_sparsity_pattern);
166  }
167 
168 
169 
170  template <>
171  void
173  reinit (const BlockSparsityPattern &block_sparsity_pattern)
174  {
175 
176  // Call the other basic reinit function, ...
177  reinit (block_sparsity_pattern.n_block_rows(),
178  block_sparsity_pattern.n_block_cols());
179 
180  // ... set the correct sizes, ...
181  this->row_block_indices = block_sparsity_pattern.get_row_indices();
182  this->column_block_indices = block_sparsity_pattern.get_column_indices();
183 
184  // ... and then assign the correct
185  // data to the blocks.
186  for (size_type r=0; r<this->n_block_rows(); ++r)
187  for (size_type c=0; c<this->n_block_cols(); ++c)
188  {
189  this->sub_objects[r][c]->reinit (block_sparsity_pattern.block(r,c));
190  }
191  }
192 
193 
194 
195  void
197  reinit (const std::vector<Epetra_Map> &parallel_partitioning,
198  const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
199  const double drop_tolerance)
200  {
201  const size_type n_block_rows = parallel_partitioning.size();
202 
203  Assert (n_block_rows == dealii_block_sparse_matrix.n_block_rows(),
205  dealii_block_sparse_matrix.n_block_rows()));
206  Assert (n_block_rows == dealii_block_sparse_matrix.n_block_cols(),
208  dealii_block_sparse_matrix.n_block_cols()));
209 
210  // Call the other basic reinit function ...
212 
213  // ... and then assign the correct
214  // data to the blocks.
215  for (size_type r=0; r<this->n_block_rows(); ++r)
216  for (size_type c=0; c<this->n_block_cols(); ++c)
217  {
218  this->sub_objects[r][c]->reinit(parallel_partitioning[r],
219  parallel_partitioning[c],
220  dealii_block_sparse_matrix.block(r,c),
221  drop_tolerance);
222  }
223 
224  collect_sizes();
225  }
226 
227 
228 
229  void
231  reinit (const ::BlockSparseMatrix<double> &dealii_block_sparse_matrix,
232  const double drop_tolerance)
233  {
234  Assert (dealii_block_sparse_matrix.n_block_rows() ==
235  dealii_block_sparse_matrix.n_block_cols(),
236  ExcDimensionMismatch (dealii_block_sparse_matrix.n_block_rows(),
237  dealii_block_sparse_matrix.n_block_cols()));
238  Assert (dealii_block_sparse_matrix.m() ==
239  dealii_block_sparse_matrix.n(),
240  ExcDimensionMismatch (dealii_block_sparse_matrix.m(),
241  dealii_block_sparse_matrix.n()));
242 
243  // produce a dummy local map and pass it
244  // off to the other function
245 #ifdef DEAL_II_WITH_MPI
246  Epetra_MpiComm trilinos_communicator (MPI_COMM_SELF);
247 #else
248  Epetra_SerialComm trilinos_communicator;
249 #endif
250 
251  std::vector<Epetra_Map> parallel_partitioning;
252  for (size_type i=0; i<dealii_block_sparse_matrix.n_block_rows(); ++i)
253  parallel_partitioning.push_back (Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(dealii_block_sparse_matrix.block(i,0).m()),
254  0,
255  trilinos_communicator));
256 
257  reinit (parallel_partitioning, dealii_block_sparse_matrix, drop_tolerance);
258  }
259 
260 
261 
262 
263 
264  void
266  {
267  // simply forward to the (non-public) function of the base class
269  }
270 
271 
272 
273  BlockSparseMatrix::size_type
275  {
276  size_type n_nonzero = 0;
277  for (size_type rows = 0; rows<this->n_block_rows(); ++rows)
278  for (size_type cols = 0; cols<this->n_block_cols(); ++cols)
279  n_nonzero += this->block(rows,cols).n_nonzero_elements();
280 
281  return n_nonzero;
282  }
283 
284 
285 
286  TrilinosScalar
288  const MPI::BlockVector &x,
289  const MPI::BlockVector &b) const
290  {
291  vmult (dst, x);
292  dst -= b;
293  dst *= -1.;
294 
295  return dst.l2_norm();
296  }
297 
298 
299 
300  // TODO: In the following we
301  // use the same code as just
302  // above six more times. Use
303  // templates.
304  TrilinosScalar
306  const BlockVector &x,
307  const BlockVector &b) const
308  {
309  vmult (dst, x);
310  dst -= b;
311  dst *= -1.;
312 
313  return dst.l2_norm();
314  }
315 
316 
317 
318  TrilinosScalar
320  const MPI::Vector &x,
321  const MPI::BlockVector &b) const
322  {
323  vmult (dst, x);
324  dst -= b;
325  dst *= -1.;
326 
327  return dst.l2_norm();
328  }
329 
330 
331 
332  TrilinosScalar
334  const Vector &x,
335  const BlockVector &b) const
336  {
337  vmult (dst, x);
338  dst -= b;
339  dst *= -1.;
340 
341  return dst.l2_norm();
342  }
343 
344 
345 
346  TrilinosScalar
348  const MPI::BlockVector &x,
349  const MPI::Vector &b) const
350  {
351  vmult (dst, x);
352  dst -= b;
353  dst *= -1.;
354 
355  return dst.l2_norm();
356  }
357 
358 
359 
360  TrilinosScalar
362  const BlockVector &x,
363  const Vector &b) const
364  {
365  vmult (dst, x);
366  dst -= b;
367  dst *= -1.;
368 
369  return dst.l2_norm();
370  }
371 
372 
373 
374  TrilinosScalar
376  const VectorBase &x,
377  const VectorBase &b) const
378  {
379  vmult (dst, x);
380  dst -= b;
381  dst *= -1.;
382 
383  return dst.l2_norm();
384  }
385 
386 
387 
388  std::vector<Epetra_Map>
390  {
391  Assert (this->n_block_cols() != 0, ExcNotInitialized());
392  Assert (this->n_block_rows() != 0, ExcNotInitialized());
393 
394  std::vector<Epetra_Map> domain_partitioner;
395  for (size_type c = 0; c < this->n_block_cols(); ++c)
396  domain_partitioner.push_back(this->sub_objects[0][c]->domain_partitioner());
397 
398  return domain_partitioner;
399  }
400 
401 
402 
403  std::vector<Epetra_Map>
405  {
406  Assert (this->n_block_cols() != 0, ExcNotInitialized());
407  Assert (this->n_block_rows() != 0, ExcNotInitialized());
408 
409  std::vector<Epetra_Map> range_partitioner;
410  for (size_type r = 0; r < this->n_block_rows(); ++r)
411  range_partitioner.push_back(this->sub_objects[r][0]->range_partitioner());
412 
413  return range_partitioner;
414  }
415 
416 
417 
418 
419 
420 
421 
422  // -------------------- explicit instantiations -----------------------
423  //
424  template void
425  BlockSparseMatrix::reinit (const ::BlockSparsityPattern &);
426  template void
427  BlockSparseMatrix::reinit (const ::BlockDynamicSparsityPattern &);
428 
429  template void
430  BlockSparseMatrix::reinit (const std::vector<Epetra_Map> &,
431  const ::BlockSparsityPattern &,
432  const bool);
433  template void
434  BlockSparseMatrix::reinit (const std::vector<Epetra_Map> &,
435  const ::BlockDynamicSparsityPattern &,
436  const bool);
437 
438  template void
439  BlockSparseMatrix::reinit (const std::vector<IndexSet> &,
440  const ::BlockDynamicSparsityPattern &,
441  const MPI_Comm &,
442  const bool);
443 
444 }
445 
446 
447 DEAL_II_NAMESPACE_CLOSE
448 
449 #endif
real_type l2_norm() const
Subscriptor & operator=(const Subscriptor &)
Definition: subscriptor.cc:154
static ::ExceptionBase & ExcNotInitialized()
const Epetra_Comm & comm_self()
Definition: utilities.cc:736
SparsityPatternType & block(const size_type row, const size_type column)
std::size_t n_nonzero_elements() const
BlockSparseMatrix & operator=(const BlockSparseMatrix &)
unsigned int n_block_cols() const
void vmult(VectorType1 &dst, const VectorType2 &src) const
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< Epetra_Map > range_partitioner() const 1
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
void reinit(const size_type n_block_rows, const size_type n_block_columns)
const BlockIndices & get_row_indices() const
BlockType & block(const unsigned int row, const unsigned int column)
Table< 2, SmartPointer< BlockType, BlockMatrixBase< SparseMatrix > > > sub_objects
std::vector< Epetra_Map > domain_partitioner() const 1
unsigned int n_block_rows() const
const BlockIndices & get_column_indices() const
TrilinosScalar residual(MPI::BlockVector &dst, const MPI::BlockVector &x, const MPI::BlockVector &b) const
real_type l2_norm() const
static ::ExceptionBase & ExcInternalError()