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