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