Reference documentation for deal.II version 9.0.0
block_matrix_array.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 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 
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>
22 
23 #include <deal.II/lac/petsc_parallel_block_vector.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 template <typename number, typename BlockVectorType>
30  :
31  row(e.row),
32  col(e.col),
33  prefix(e.prefix),
35  matrix(e.matrix)
36 {
37  Entry &e2 = const_cast<Entry &>(e);
38  e2.matrix = nullptr;
39 }
40 
41 
42 
43 template <typename number, typename BlockVectorType>
45 {
46  if (matrix)
47  delete matrix;
48 }
49 
50 
51 template <typename number, typename BlockVectorType>
53  : block_rows (0),
54  block_cols (0)
55 {}
56 
57 
58 
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)
65 {}
66 
67 
68 template <typename number, typename BlockVectorType>
69 void
71 (const unsigned int n_block_rows,
72  const unsigned int n_block_cols)
73 {
74  block_rows = n_block_rows;
75  block_cols = n_block_cols;
76 }
77 
78 
79 
80 template <typename number, typename BlockVectorType>
81 void
83 (const unsigned int n_block_rows,
84  const unsigned int n_block_cols)
85 {
86  clear();
87  block_rows = n_block_rows;
88  block_cols = n_block_cols;
89 }
90 
91 
92 
93 template <typename number, typename BlockVectorType>
94 void
96 {
97  entries.clear();
98 }
99 
100 
101 template <typename number, typename BlockVectorType>
102 void
104  const BlockVectorType &src) const
105 {
107  Assert (dst.n_blocks() == block_rows,
108  ExcDimensionMismatch(dst.n_blocks(), block_rows));
109  Assert (src.n_blocks() == block_cols,
110  ExcDimensionMismatch(src.n_blocks(), block_cols));
111 
113  typename BlockVectorType::BlockType &aux = *p_aux;
114 
115  typename std::vector<Entry>::const_iterator m = entries.begin();
116  typename std::vector<Entry>::const_iterator end = entries.end();
117 
118  for (; m != end ; ++m)
119  {
120  aux.reinit(dst.block(m->row));
121  if (m->transpose)
122  m->matrix->Tvmult(aux, src.block(m->col));
123  else
124  m->matrix->vmult(aux, src.block(m->col));
125  dst.block(m->row).add (m->prefix, aux);
126  }
127 }
128 
129 
130 
131 
132 template <typename number, typename BlockVectorType>
133 void
135  const BlockVectorType &src) const
136 {
137  dst = 0.;
138  vmult_add (dst, src);
139 }
140 
141 
142 
143 
144 template <typename number, typename BlockVectorType>
145 void
147  const BlockVectorType &src) const
148 {
150  Assert (dst.n_blocks() == block_cols,
151  ExcDimensionMismatch(dst.n_blocks(), block_cols));
152  Assert (src.n_blocks() == block_rows,
153  ExcDimensionMismatch(src.n_blocks(), block_rows));
154 
155  typename std::vector<Entry>::const_iterator m = entries.begin();
156  typename std::vector<Entry>::const_iterator end = entries.end();
157 
159  typename BlockVectorType::BlockType &aux = *p_aux;
160 
161  for (; m != end ; ++m)
162  {
163  aux.reinit(dst.block(m->col));
164  if (m->transpose)
165  m->matrix->vmult(aux, src.block(m->row));
166  else
167  m->matrix->Tvmult(aux, src.block(m->row));
168  dst.block(m->col).add (m->prefix, aux);
169  }
170 }
171 
172 
173 
174 template <typename number, typename BlockVectorType>
175 void
177  const BlockVectorType &src) const
178 {
179  dst = 0.;
180  Tvmult_add (dst, src);
181 }
182 
183 
184 
185 
186 template <typename number, typename BlockVectorType>
187 number
189 (const BlockVectorType &u,
190  const BlockVectorType &v) const
191 {
193  Assert (u.n_blocks() == block_rows,
194  ExcDimensionMismatch(u.n_blocks(), block_rows));
195  Assert (v.n_blocks() == block_cols,
196  ExcDimensionMismatch(v.n_blocks(), block_cols));
197 
199  typename BlockVectorType::BlockType &aux = *p_aux;
200 
201  typename std::vector<Entry>::const_iterator m;
202  typename std::vector<Entry>::const_iterator end = entries.end();
203 
204  number result = 0.;
205 
206  for (unsigned int i=0; i<block_rows; ++i)
207  {
208  aux.reinit(u.block(i));
209  for (m = entries.begin(); m != end ; ++m)
210  {
211  if (m->row != i)
212  continue;
213  if (m->transpose)
214  m->matrix->Tvmult_add(aux, v.block(m->col));
215  else
216  m->matrix->vmult(aux, v.block(m->col));
217  }
218  result += u.block(i)*aux;
219  }
220 
221  return result;
222 }
223 
224 
225 
226 template <typename number, typename BlockVectorType>
227 number
229 (const BlockVectorType &u) const
230 {
231  return matrix_scalar_product(u,u);
232 }
233 
234 
235 
236 template <typename number, typename BlockVectorType>
237 unsigned int
239 {
240  return block_rows;
241 }
242 
243 
244 
245 template <typename number, typename BlockVectorType>
246 unsigned int
248 {
249  return block_cols;
250 }
251 
252 
253 
254 //---------------------------------------------------------------------------
255 
256 template <typename number, typename BlockVectorType>
258  : BlockMatrixArray<number,BlockVectorType> (),
259  backward(false)
260 {}
261 
262 
263 template <typename number, typename BlockVectorType>
265 (const unsigned int block_rows)
266  :
267  BlockMatrixArray<number,BlockVectorType> (block_rows, block_rows),
268  backward(false)
269 {}
270 
271 
272 template <typename number, typename BlockVectorType>
273 void
275 (const unsigned int n)
276 {
278 }
279 
280 
281 template <typename number, typename BlockVectorType>
282 void
284 (BlockVectorType &dst,
285  size_type row_num) const
286 {
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>
293  diagonals;
294 
296  typename BlockVectorType::BlockType &aux = *p_aux;
297 
298  aux.reinit(dst.block(row_num), true);
299 
300  // Loop over all entries, since
301  // they are not ordered by rows.
302  for (; m != end ; ++m)
303  {
304  const size_type i=m->row;
305  // Ignore everything not in
306  // this row
307  if (i != row_num)
308  continue;
309  const size_type j=m->col;
310  // Only use the lower (upper)
311  // triangle for forward
312  // (backward) substitution
313  if (((j > i) && !backward) || ((j < i) && backward))
314  continue;
315  if (j == i)
316  {
317  diagonals.push_back(m);
318  }
319  else
320  {
321  if (m->transpose)
322  m->matrix->Tvmult(aux, dst.block(j));
323  else
324  m->matrix->vmult(aux, dst.block(j));
325  dst.block(i).add (-1.0 * m->prefix, aux);
326  }
327  }
328  Assert (diagonals.size() != 0, ExcNoDiagonal(row_num));
329 
330  // Inverting the diagonal block is
331  // simple, if there is only one
332  // matrix
333  if (diagonals.size() == 1)
334  {
335  if (diagonals[0]->transpose)
336  diagonals[0]->matrix->Tvmult(aux, dst.block(row_num));
337  else
338  diagonals[0]->matrix->vmult(aux, dst.block(row_num));
339  dst.block(row_num).equ (diagonals[0]->prefix, aux);
340  }
341  else
342  {
343  aux = 0.;
344  for (size_type i=0; i<diagonals.size(); ++i)
345  {
346  m = diagonals[i];
347  // First, divide by the current
348  // factor, such that we can
349  // multiply by it later.
350  aux /= m->prefix;
351  if (m->transpose)
352  m->matrix->Tvmult_add(aux, dst.block(row_num));
353  else
354  m->matrix->vmult_add(aux, dst.block(row_num));
355  aux *= m->prefix;
356  }
357  dst.block(row_num) = aux;
358  }
359 }
360 
361 
362 
363 template <typename number, typename BlockVectorType>
364 void
366 (BlockVectorType &dst,
367  const BlockVectorType &src) const
368 {
369  Assert (dst.n_blocks() == n_block_rows(),
370  ExcDimensionMismatch(dst.n_blocks(), n_block_rows()));
371  Assert (src.n_blocks() == n_block_cols(),
372  ExcDimensionMismatch(src.n_blocks(), n_block_cols()));
373 
374  BlockVectorType aux;
375  aux.reinit(dst);
376  vmult(aux, src);
377  dst += aux;
378 }
379 
380 
381 
382 template <typename number, typename BlockVectorType>
383 void
385 (BlockVectorType &dst,
386  const BlockVectorType &src) const
387 {
388  Assert (dst.n_blocks() == n_block_rows(),
389  ExcDimensionMismatch(dst.n_blocks(), n_block_rows()));
390  Assert (src.n_blocks() == n_block_cols(),
391  ExcDimensionMismatch(src.n_blocks(), n_block_cols()));
392 
393  dst.equ(1., src);
394 
395  if (backward)
396  {
397  for (unsigned int i=n_block_rows(); i>0;)
398  do_row(dst, --i);
399  }
400  else
401  {
402  for (unsigned int i=0; i<n_block_rows(); ++i)
403  do_row(dst, i);
404  }
405 
406 }
407 
408 template <typename number, typename BlockVectorType>
409 void
411 (BlockVectorType &,
412  const BlockVectorType &) const
413 {
414  Assert (false, ExcNotImplemented());
415 }
416 
417 
418 template <typename number, typename BlockVectorType>
419 void
421 (BlockVectorType &,
422  const BlockVectorType &) const
423 {
424  Assert (false, ExcNotImplemented());
425 }
426 
427 template class BlockMatrixArray<float>;
428 template class BlockMatrixArray<double>;
429 template class BlockTrianglePrecondition<float>;
430 template class BlockTrianglePrecondition<double>;
431 
432 #ifdef DEAL_II_WITH_TRILINOS
437 #endif
438 
439 #ifdef DEAL_II_WITH_PETSC
440 #ifdef PETSC_USE_COMPLEX
445 #else
450 #endif
451 #endif
452 
453 DEAL_II_NAMESPACE_CLOSE
unsigned int block_cols
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)
Definition: exceptions.h:1142
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
unsigned int block_rows
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)
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
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)