Reference documentation for deal.II version 8.5.1
matrix_tools.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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/base/function.h>
17 #include <deal.II/base/quadrature.h>
18 #include <deal.II/base/work_stream.h>
19 #include <deal.II/base/geometry_info.h>
20 #include <deal.II/base/quadrature.h>
21 #include <deal.II/dofs/dof_handler.h>
22 #include <deal.II/dofs/dof_accessor.h>
23 #include <deal.II/dofs/dof_tools.h>
24 #include <deal.II/fe/fe.h>
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/fe/mapping_q1.h>
27 #include <deal.II/grid/tria_iterator.h>
28 #include <deal.II/hp/fe_values.h>
29 #include <deal.II/hp/mapping_collection.h>
30 #include <deal.II/numerics/matrix_tools.h>
31 #include <deal.II/lac/vector.h>
32 #include <deal.II/lac/block_vector.h>
33 #include <deal.II/lac/full_matrix.h>
34 #include <deal.II/lac/sparse_matrix.h>
35 #include <deal.II/lac/block_sparse_matrix.h>
36 
37 #ifdef DEAL_II_WITH_PETSC
38 # include <deal.II/lac/petsc_parallel_sparse_matrix.h>
39 # include <deal.II/lac/petsc_sparse_matrix.h>
40 # include <deal.II/lac/petsc_parallel_vector.h>
41 # include <deal.II/lac/petsc_vector.h>
42 # include <deal.II/lac/petsc_parallel_block_sparse_matrix.h>
43 #endif
44 
45 #ifdef DEAL_II_WITH_TRILINOS
46 # include <deal.II/lac/trilinos_sparse_matrix.h>
47 # include <deal.II/lac/trilinos_vector.h>
48 # include <deal.II/lac/trilinos_block_sparse_matrix.h>
49 # include <deal.II/lac/trilinos_block_vector.h>
50 #endif
51 
52 #include <algorithm>
53 
54 
55 #include <algorithm>
56 #include <set>
57 #include <cmath>
58 
59 
60 DEAL_II_NAMESPACE_OPEN
61 
62 
63 
64 
65 namespace MatrixTools
66 {
67  namespace
68  {
69  template <typename Iterator>
70  bool column_less_than(const typename Iterator::value_type p,
71  const unsigned int column)
72  {
73  return (p.column() < column);
74  }
75  }
76 
77 //TODO:[WB] I don't think that the optimized storage of diagonals is needed (GK)
78  template <typename number>
79  void
80  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
81  SparseMatrix<number> &matrix,
82  Vector<number> &solution,
83  Vector<number> &right_hand_side,
84  const bool eliminate_columns)
85  {
86  Assert (matrix.n() == right_hand_side.size(),
87  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
88  Assert (matrix.n() == solution.size(),
89  ExcDimensionMismatch(matrix.n(), solution.size()));
90  Assert (matrix.n() == matrix.m(),
91  ExcDimensionMismatch(matrix.n(), matrix.m()));
92 
93  // if no boundary values are to be applied
94  // simply return
95  if (boundary_values.size() == 0)
96  return;
97 
98 
99  const types::global_dof_index n_dofs = matrix.m();
100 
101  // if a diagonal entry is zero
102  // later, then we use another
103  // number instead. take it to be
104  // the first nonzero diagonal
105  // element of the matrix, or 1 if
106  // there is no such thing
107  number first_nonzero_diagonal_entry = 1;
108  for (unsigned int i=0; i<n_dofs; ++i)
109  if (matrix.diag_element(i) != number())
110  {
111  first_nonzero_diagonal_entry = matrix.diag_element(i);
112  break;
113  }
114 
115 
116  typename std::map<types::global_dof_index,number>::const_iterator dof = boundary_values.begin(),
117  endd = boundary_values.end();
118  for (; dof != endd; ++dof)
119  {
120  Assert (dof->first < n_dofs, ExcInternalError());
121 
122  const types::global_dof_index dof_number = dof->first;
123  // for each boundary dof:
124 
125  // set entries of this line to zero except for the diagonal
126  // entry
127  for (typename SparseMatrix<number>::iterator
128  p = matrix.begin(dof_number);
129  p != matrix.end(dof_number); ++p)
130  if (p->column() != dof_number)
131  p->value() = 0.;
132 
133  // set right hand side to
134  // wanted value: if main diagonal
135  // entry nonzero, don't touch it
136  // and scale rhs accordingly. If
137  // zero, take the first main
138  // diagonal entry we can find, or
139  // one if no nonzero main diagonal
140  // element exists. Normally, however,
141  // the main diagonal entry should
142  // not be zero.
143  //
144  // store the new rhs entry to make
145  // the gauss step more efficient
146  number new_rhs;
147  if (matrix.diag_element(dof_number) != number())
148  {
149  new_rhs = dof->second * matrix.diag_element(dof_number);
150  right_hand_side(dof_number) = new_rhs;
151  }
152  else
153  {
154  matrix.set (dof_number, dof_number,
155  first_nonzero_diagonal_entry);
156  new_rhs = dof->second * first_nonzero_diagonal_entry;
157  right_hand_side(dof_number) = new_rhs;
158  }
159 
160 
161  // if the user wants to have
162  // the symmetry of the matrix
163  // preserved, and if the
164  // sparsity pattern is
165  // symmetric, then do a Gauss
166  // elimination step with the
167  // present row
168  if (eliminate_columns)
169  {
170  // store the only nonzero entry
171  // of this line for the Gauss
172  // elimination step
173  const number diagonal_entry = matrix.diag_element(dof_number);
174 
175  // we have to loop over all rows of the matrix which have
176  // a nonzero entry in the column which we work in
177  // presently. if the sparsity pattern is symmetric, then
178  // we can get the positions of these rows cheaply by
179  // looking at the nonzero column numbers of the present
180  // row. we need not look at the first entry of each row,
181  // since that is the diagonal element and thus the present
182  // row
183  for (typename SparseMatrix<number>::iterator
184  q = matrix.begin(dof_number)+1;
185  q != matrix.end(dof_number); ++q)
186  {
187  const types::global_dof_index row = q->column();
188 
189  // find the position of
190  // element
191  // (row,dof_number)
192  bool (*comp)(const typename SparseMatrix<number>::iterator::value_type p,
193  const unsigned int column)
194  = &column_less_than<typename SparseMatrix<number>::iterator>;
195  const typename SparseMatrix<number>::iterator
196  p = Utilities::lower_bound(matrix.begin(row)+1,
197  matrix.end(row),
198  dof_number,
199  comp);
200 
201  // check whether this line has an entry in the
202  // regarding column (check for ==dof_number and !=
203  // next_row, since if row==dof_number-1, *p is a
204  // past-the-end pointer but points to dof_number
205  // anyway...)
206  //
207  // there should be such an entry! we know this because
208  // we have assumed that the sparsity pattern is
209  // symmetric and we only walk over those rows for
210  // which the current row has a column entry
211  Assert ((p != matrix.end(row))
212  &&
213  (p->column() == dof_number),
214  ExcMessage("This function is trying to access an element of the "
215  "matrix that doesn't seem to exist. Are you using a "
216  "nonsymmetric sparsity pattern? If so, you are not "
217  "allowed to set the eliminate_column argument of this "
218  "function, see the documentation."));
219 
220  // correct right hand side
221  right_hand_side(row) -= static_cast<number>(p->value()) /
222  diagonal_entry * new_rhs;
223 
224  // set matrix entry to zero
225  p->value() = 0.;
226  }
227  }
228 
229  // preset solution vector
230  solution(dof_number) = dof->second;
231  }
232  }
233 
234 
235 
236  template <typename number>
237  void
238  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
240  BlockVector<number> &solution,
241  BlockVector<number> &right_hand_side,
242  const bool eliminate_columns)
243  {
244  const unsigned int blocks = matrix.n_block_rows();
245 
246  Assert (matrix.n() == right_hand_side.size(),
247  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
248  Assert (matrix.n() == solution.size(),
249  ExcDimensionMismatch(matrix.n(), solution.size()));
250  Assert (matrix.n_block_rows() == matrix.n_block_cols(),
251  ExcNotQuadratic());
252  Assert (matrix.get_sparsity_pattern().get_row_indices() ==
253  matrix.get_sparsity_pattern().get_column_indices(),
254  ExcNotQuadratic());
255  Assert (matrix.get_sparsity_pattern().get_column_indices() ==
256  solution.get_block_indices (),
257  ExcBlocksDontMatch ());
258  Assert (matrix.get_sparsity_pattern().get_row_indices() ==
259  right_hand_side.get_block_indices (),
260  ExcBlocksDontMatch ());
261 
262  // if no boundary values are to be applied
263  // simply return
264  if (boundary_values.size() == 0)
265  return;
266 
267 
268  const types::global_dof_index n_dofs = matrix.m();
269 
270  // if a diagonal entry is zero
271  // later, then we use another
272  // number instead. take it to be
273  // the first nonzero diagonal
274  // element of the matrix, or 1 if
275  // there is no such thing
276  number first_nonzero_diagonal_entry = 0;
277  for (unsigned int diag_block=0; diag_block<blocks; ++diag_block)
278  {
279  for (unsigned int i=0; i<matrix.block(diag_block,diag_block).n(); ++i)
280  if (matrix.block(diag_block,diag_block).diag_element(i) != 0)
281  {
282  first_nonzero_diagonal_entry
283  = matrix.block(diag_block,diag_block).diag_element(i);
284  break;
285  }
286  // check whether we have found
287  // something in the present
288  // block
289  if (first_nonzero_diagonal_entry != 0)
290  break;
291  }
292  // nothing found on all diagonal
293  // blocks? if so, use 1.0 instead
294  if (first_nonzero_diagonal_entry == 0)
295  first_nonzero_diagonal_entry = 1;
296 
297 
298  typename std::map<types::global_dof_index,number>::const_iterator dof = boundary_values.begin(),
299  endd = boundary_values.end();
300  const BlockSparsityPattern &
301  sparsity_pattern = matrix.get_sparsity_pattern();
302 
303  // pointer to the mapping between
304  // global and block indices. since
305  // the row and column mappings are
306  // equal, store a pointer on only
307  // one of them
308  const BlockIndices &
309  index_mapping = sparsity_pattern.get_column_indices();
310 
311  // now loop over all boundary dofs
312  for (; dof != endd; ++dof)
313  {
314  Assert (dof->first < n_dofs, ExcInternalError());
315  (void)n_dofs;
316 
317  // get global index and index
318  // in the block in which this
319  // dof is located
320  const types::global_dof_index dof_number = dof->first;
321  const std::pair<unsigned int,types::global_dof_index>
322  block_index = index_mapping.global_to_local (dof_number);
323 
324  // for each boundary dof:
325 
326  // set entries of this line
327  // to zero except for the diagonal
328  // entry. Note that the diagonal
329  // entry is always the first one
330  // in a row for square matrices
331  for (unsigned int block_col=0; block_col<blocks; ++block_col)
332  for (typename SparseMatrix<number>::iterator
333  p = (block_col == block_index.first ?
334  matrix.block(block_index.first,block_col).begin(block_index.second) + 1 :
335  matrix.block(block_index.first,block_col).begin(block_index.second));
336  p != matrix.block(block_index.first,block_col).end(block_index.second);
337  ++p)
338  p->value() = 0;
339 
340  // set right hand side to
341  // wanted value: if main diagonal
342  // entry nonzero, don't touch it
343  // and scale rhs accordingly. If
344  // zero, take the first main
345  // diagonal entry we can find, or
346  // one if no nonzero main diagonal
347  // element exists. Normally, however,
348  // the main diagonal entry should
349  // not be zero.
350  //
351  // store the new rhs entry to make
352  // the gauss step more efficient
353  number new_rhs;
354  if (matrix.block(block_index.first, block_index.first)
355  .diag_element(block_index.second) != 0.0)
356  new_rhs = dof->second *
357  matrix.block(block_index.first, block_index.first)
358  .diag_element(block_index.second);
359  else
360  {
361  matrix.block(block_index.first, block_index.first)
362  .diag_element(block_index.second)
363  = first_nonzero_diagonal_entry;
364  new_rhs = dof->second * first_nonzero_diagonal_entry;
365  }
366  right_hand_side.block(block_index.first)(block_index.second)
367  = new_rhs;
368 
369 
370  // if the user wants to have
371  // the symmetry of the matrix
372  // preserved, and if the
373  // sparsity pattern is
374  // symmetric, then do a Gauss
375  // elimination step with the
376  // present row. this is a
377  // little more complicated for
378  // block matrices.
379  if (eliminate_columns)
380  {
381  // store the only nonzero entry
382  // of this line for the Gauss
383  // elimination step
384  const number diagonal_entry
385  = matrix.block(block_index.first,block_index.first)
386  .diag_element(block_index.second);
387 
388  // we have to loop over all
389  // rows of the matrix which
390  // have a nonzero entry in
391  // the column which we work
392  // in presently. if the
393  // sparsity pattern is
394  // symmetric, then we can
395  // get the positions of
396  // these rows cheaply by
397  // looking at the nonzero
398  // column numbers of the
399  // present row.
400  //
401  // note that if we check
402  // whether row @p{row} in
403  // block (r,c) is non-zero,
404  // then we have to check
405  // for the existence of
406  // column @p{row} in block
407  // (c,r), i.e. of the
408  // transpose block
409  for (unsigned int block_row=0; block_row<blocks; ++block_row)
410  {
411  // get pointers to the sparsity patterns of this block and of
412  // the transpose one
413  const SparsityPattern &this_sparsity
414  = sparsity_pattern.block (block_row, block_index.first);
415 
416  SparseMatrix<number> &this_matrix
417  = matrix.block(block_row, block_index.first);
418  SparseMatrix<number> &transpose_matrix
419  = matrix.block(block_index.first, block_row);
420 
421  // traverse the row of the transpose block to find the
422  // interesting rows in the present block. don't use the
423  // diagonal element of the diagonal block
424  for (typename SparseMatrix<number>::iterator
425  q = (block_index.first == block_row ?
426  transpose_matrix.begin(block_index.second)+1 :
427  transpose_matrix.begin(block_index.second));
428  q != transpose_matrix.end(block_index.second);
429  ++q)
430  {
431  // get the number of the column in this row in which a
432  // nonzero entry is. this is also the row of the transpose
433  // block which has an entry in the interesting row
434  const types::global_dof_index row = q->column();
435 
436  // find the position of element (row,dof_number) in this
437  // block (not in the transpose one). note that we have to
438  // take care of special cases with square sub-matrices
439  bool (*comp)(typename SparseMatrix<number>::iterator::value_type p,
440  const unsigned int column)
441  = &column_less_than<typename SparseMatrix<number>::iterator>;
442 
443  typename SparseMatrix<number>::iterator p = this_matrix.end();
444 
445  if (this_sparsity.n_rows() == this_sparsity.n_cols())
446  {
447  if (this_matrix.begin(row)->column()
448  ==
449  block_index.second)
450  p = this_matrix.begin(row);
451  else
452  p = Utilities::lower_bound(this_matrix.begin(row)+1,
453  this_matrix.end(row),
454  block_index.second,
455  comp);
456  }
457  else
458  p = Utilities::lower_bound(this_matrix.begin(row),
459  this_matrix.end(row),
460  block_index.second,
461  comp);
462 
463  // check whether this line has an entry in the
464  // regarding column (check for ==dof_number and !=
465  // next_row, since if row==dof_number-1, *p is a
466  // past-the-end pointer but points to dof_number
467  // anyway...)
468  //
469  // there should be such an entry! we know this because
470  // we have assumed that the sparsity pattern is
471  // symmetric and we only walk over those rows for
472  // which the current row has a column entry
473  Assert ((p->column() == block_index.second) &&
474  (p != this_matrix.end(row)),
475  ExcInternalError());
476 
477  // correct right hand side
478  right_hand_side.block(block_row)(row)
479  -= p->value() /
480  diagonal_entry * new_rhs;
481 
482  // set matrix entry to zero
483  p->value() = 0.;
484  }
485  }
486  }
487 
488  // preset solution vector
489  solution.block(block_index.first)(block_index.second) = dof->second;
490  }
491  }
492 
493 
494 
495 
496 
497  template <typename number>
498  void
499  local_apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
500  const std::vector<types::global_dof_index> &local_dof_indices,
501  FullMatrix<number> &local_matrix,
502  Vector<number> &local_rhs,
503  const bool eliminate_columns)
504  {
505  Assert (local_dof_indices.size() == local_matrix.m(),
506  ExcDimensionMismatch(local_dof_indices.size(),
507  local_matrix.m()));
508  Assert (local_dof_indices.size() == local_matrix.n(),
509  ExcDimensionMismatch(local_dof_indices.size(),
510  local_matrix.n()));
511  Assert (local_dof_indices.size() == local_rhs.size(),
512  ExcDimensionMismatch(local_dof_indices.size(),
513  local_rhs.size()));
514 
515  // if there is nothing to do, then exit
516  // right away
517  if (boundary_values.size() == 0)
518  return;
519 
520  // otherwise traverse all the dofs used in
521  // the local matrices and vectors and see
522  // what's there to do
523 
524  // if we need to treat an entry, then we
525  // set the diagonal entry to its absolute
526  // value. if it is zero, we used to set it
527  // to one, which is a really terrible
528  // choice that can lead to hours of
529  // searching for bugs in programs (I
530  // experienced this :-( ) if the matrix
531  // entries are otherwise very large. this
532  // is so since iterative solvers would
533  // simply not correct boundary nodes for
534  // their correct values since the residual
535  // contributions of their rows of the
536  // linear system is almost zero if the
537  // diagonal entry is one. thus, set it to
538  // the average absolute value of the
539  // nonzero diagonal elements.
540  //
541  // we only compute this value lazily the
542  // first time we need it.
543  number average_diagonal = 0;
544  const unsigned int n_local_dofs = local_dof_indices.size();
545  for (unsigned int i=0; i<n_local_dofs; ++i)
546  {
547  const typename std::map<types::global_dof_index, number>::const_iterator
548  boundary_value = boundary_values.find (local_dof_indices[i]);
549  if (boundary_value != boundary_values.end())
550  {
551  // remove this row, except for the
552  // diagonal element
553  for (unsigned int j=0; j<n_local_dofs; ++j)
554  if (i != j)
555  local_matrix(i,j) = 0;
556 
557  // replace diagonal entry by its
558  // absolute value to make sure that
559  // everything remains positive, or
560  // by the average diagonal value if
561  // zero
562  if (local_matrix(i,i) == 0.)
563  {
564  // if average diagonal hasn't
565  // yet been computed, do so now
566  if (average_diagonal == 0.)
567  {
568  unsigned int nonzero_diagonals = 0;
569  for (unsigned int k=0; k<n_local_dofs; ++k)
570  if (local_matrix(k,k) != 0.)
571  {
572  average_diagonal += std::fabs(local_matrix(k,k));
573  ++nonzero_diagonals;
574  }
575  if (nonzero_diagonals != 0)
576  average_diagonal /= nonzero_diagonals;
577  else
578  average_diagonal = 0;
579  }
580 
581  // only if all diagonal entries
582  // are zero, then resort to the
583  // last measure: choose one
584  if (average_diagonal == 0.)
585  average_diagonal = 1.;
586 
587  local_matrix(i,i) = average_diagonal;
588  }
589  else
590  local_matrix(i,i) = std::fabs(local_matrix(i,i));
591 
592  // and replace rhs entry by correct
593  // value
594  local_rhs(i) = local_matrix(i,i) * boundary_value->second;
595 
596  // finally do the elimination step
597  // if requested
598  if (eliminate_columns == true)
599  {
600  for (unsigned int row=0; row<n_local_dofs; ++row)
601  if (row != i)
602  {
603  local_rhs(row) -= local_matrix(row,i) * boundary_value->second;
604  local_matrix(row,i) = 0;
605  }
606  }
607  }
608  }
609  }
610 }
611 
612 
613 
614 // explicit instantiations
615 #include "matrix_tools.inst"
616 
617 
618 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:618
size_type m() const
void local_apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, const std::vector< types::global_dof_index > &local_dof_indices, FullMatrix< number > &local_matrix, Vector< number > &local_rhs, const bool eliminate_columns)
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Definition: matrix_tools.cc:80
const_iterator end() const
size_type n_cols() const
SparsityPatternType & block(const size_type row, const size_type column)
size_type n_rows() const
static ::ExceptionBase & ExcMessage(std::string arg1)
std::size_t size() const
size_type n() const
unsigned int global_dof_index
Definition: types.h:88
const BlockIndices & get_block_indices() const
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t size() const
const Accessor< number, Constness > & value_type
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcBlocksDontMatch()
const_iterator begin() const
BlockType & block(const unsigned int i)
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
const BlockIndices & get_column_indices() const
static ::ExceptionBase & ExcInternalError()