Reference documentation for deal.II version 9.0.0
matrix_tools_once.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_matrix_base.h>
39 # include <deal.II/lac/petsc_parallel_block_vector.h>
40 # include <deal.II/lac/petsc_parallel_block_sparse_matrix.h>
41 # include <deal.II/lac/petsc_vector_base.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 #include <cmath>
53 
54 
55 DEAL_II_NAMESPACE_OPEN
56 
57 namespace MatrixTools
58 {
59 
60 #ifdef DEAL_II_WITH_PETSC
61  void
63  (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
65  PETScWrappers::VectorBase &solution,
66  PETScWrappers::VectorBase &right_hand_side,
67  const bool eliminate_columns)
68  {
69  (void)eliminate_columns;
70  Assert (eliminate_columns == false, ExcNotImplemented());
71 
72  Assert (matrix.n() == right_hand_side.size(),
73  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
74  Assert (matrix.n() == solution.size(),
75  ExcDimensionMismatch(matrix.n(), solution.size()));
76 
77  // if no boundary values are to be applied, then
78  // jump straight to the compress() calls that we still have
79  // to perform because they are collective operations
80  if (boundary_values.size() > 0)
81  {
82  const std::pair<types::global_dof_index, types::global_dof_index> local_range
83  = matrix.local_range();
84  Assert (local_range == right_hand_side.local_range(),
86  Assert (local_range == solution.local_range(),
88 
89  // determine the first nonzero diagonal
90  // entry from within the part of the
91  // matrix that we can see. if we can't
92  // find such an entry, take one
93  PetscScalar average_nonzero_diagonal_entry = 1;
94  for (types::global_dof_index i=local_range.first; i<local_range.second; ++i)
95  if (matrix.diag_element(i) != PetscScalar ())
96  {
97  average_nonzero_diagonal_entry = std::abs(matrix.diag_element(i));
98  break;
99  }
100 
101  // figure out which rows of the matrix we
102  // have to eliminate on this processor
103  std::vector<types::global_dof_index> constrained_rows;
104  for (std::map<types::global_dof_index,PetscScalar>::const_iterator
105  dof = boundary_values.begin();
106  dof != boundary_values.end();
107  ++dof)
108  if ((dof->first >= local_range.first) &&
109  (dof->first < local_range.second))
110  constrained_rows.push_back (dof->first);
111 
112  // then eliminate these rows and set
113  // their diagonal entry to what we have
114  // determined above. note that for petsc
115  // matrices interleaving read with write
116  // operations is very expensive. thus, we
117  // here always replace the diagonal
118  // element, rather than first checking
119  // whether it is nonzero and in that case
120  // preserving it. this is different from
121  // the case of deal.II sparse matrices
122  // treated in the other functions.
123  matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
124 
125  std::vector<types::global_dof_index> indices;
126  std::vector<PetscScalar> solution_values;
127  for (std::map<types::global_dof_index,PetscScalar>::const_iterator
128  dof = boundary_values.begin();
129  dof != boundary_values.end();
130  ++dof)
131  if ((dof->first >= local_range.first) &&
132  (dof->first < local_range.second))
133  {
134  indices.push_back (dof->first);
135  solution_values.push_back (dof->second);
136  }
137  solution.set (indices, solution_values);
138 
139  // now also set appropriate values for
140  // the rhs
141  for (unsigned int i=0; i<solution_values.size(); ++i)
142  solution_values[i] *= average_nonzero_diagonal_entry;
143 
144  right_hand_side.set (indices, solution_values);
145  }
146  else
147  {
148  // clear_rows() is a collective operation so we still have to call
149  // it:
150  std::vector<types::global_dof_index> constrained_rows;
151  matrix.clear_rows (constrained_rows, 1.);
152  }
153 
154  // clean up
156  right_hand_side.compress (VectorOperation::insert);
157  }
158 
159 
160  void
161  apply_boundary_values (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
164  PETScWrappers::MPI::BlockVector &right_hand_side,
165  const bool eliminate_columns)
166  {
167  Assert (matrix.n() == right_hand_side.size(),
168  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
169  Assert (matrix.n() == solution.size(),
170  ExcDimensionMismatch(matrix.n(), solution.size()));
171  Assert (matrix.n_block_rows() == matrix.n_block_cols(),
172  ExcNotQuadratic());
173 
174  const unsigned int n_blocks = matrix.n_block_rows();
175 
176  // We need to find the subdivision
177  // into blocks for the boundary values.
178  // To this end, generate a vector of
179  // maps with the respective indices.
180  std::vector<std::map<::types::global_dof_index,PetscScalar> > block_boundary_values(n_blocks);
181  {
182  int block = 0;
183  ::types::global_dof_index offset = 0;
184  for (std::map<types::global_dof_index,PetscScalar>::const_iterator
185  dof = boundary_values.begin();
186  dof != boundary_values.end();
187  ++dof)
188  {
189  if (dof->first >= matrix.block(block,0).m() + offset)
190  {
191  offset += matrix.block(block,0).m();
192  block++;
193  }
194  const types::global_dof_index index = dof->first - offset;
195  block_boundary_values[block].insert(std::pair<types::global_dof_index, PetscScalar> (index,dof->second));
196  }
197  }
198 
199  // Now call the non-block variants on
200  // the diagonal subblocks and the
201  // solution/rhs.
202  for (unsigned int block=0; block<n_blocks; ++block)
203  apply_boundary_values(block_boundary_values[block],
204  matrix.block(block,block),
205  solution.block(block),
206  right_hand_side.block(block),
207  eliminate_columns);
208 
209  // Finally, we need to do something
210  // about the off-diagonal matrices. This
211  // is luckily not difficult. Just clear
212  // the whole row.
213  for (unsigned int block_m=0; block_m<n_blocks; ++block_m)
214  {
215  const std::pair<types::global_dof_index, types::global_dof_index> local_range
216  = matrix.block(block_m,0).local_range();
217 
218  std::vector<types::global_dof_index> constrained_rows;
219  for (std::map<types::global_dof_index,PetscScalar>::const_iterator
220  dof = block_boundary_values[block_m].begin();
221  dof != block_boundary_values[block_m].end();
222  ++dof)
223  if ((dof->first >= local_range.first) &&
224  (dof->first < local_range.second))
225  constrained_rows.push_back (dof->first);
226 
227  for (unsigned int block_n=0; block_n<n_blocks; ++block_n)
228  if (block_m != block_n)
229  matrix.block(block_m,block_n).clear_rows(constrained_rows);
230  }
231  }
232 
233 #endif
234 
235 
236 
237 #ifdef DEAL_II_WITH_TRILINOS
238 
239  namespace internal
240  {
241  namespace TrilinosWrappers
242  {
243  template <typename TrilinosMatrix, typename TrilinosVector>
244  void
245  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
246  TrilinosMatrix &matrix,
247  TrilinosVector &solution,
248  TrilinosVector &right_hand_side,
249  const bool eliminate_columns)
250  {
251  Assert (eliminate_columns == false, ExcNotImplemented());
252  (void)eliminate_columns;
253 
254  Assert (matrix.n() == right_hand_side.size(),
255  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
256  Assert (matrix.n() == solution.size(),
257  ExcDimensionMismatch(matrix.m(), solution.size()));
258 
259  // if no boundary values are to be applied, then
260  // jump straight to the compress() calls that we still have
261  // to perform because they are collective operations
262  if (boundary_values.size() > 0)
263  {
264  const std::pair<types::global_dof_index, types::global_dof_index> local_range
265  = matrix.local_range();
266  Assert (local_range == right_hand_side.local_range(),
267  ExcInternalError());
268  Assert (local_range == solution.local_range(),
269  ExcInternalError());
270 
271  // determine the first nonzero diagonal
272  // entry from within the part of the
273  // matrix that we can see. if we can't
274  // find such an entry, take one
275  TrilinosScalar average_nonzero_diagonal_entry = 1;
276  for (types::global_dof_index i=local_range.first; i<local_range.second; ++i)
277  if (matrix.diag_element(i) != 0)
278  {
279  average_nonzero_diagonal_entry = std::fabs(matrix.diag_element(i));
280  break;
281  }
282 
283  // figure out which rows of the matrix we
284  // have to eliminate on this processor
285  std::vector<types::global_dof_index> constrained_rows;
286  for (std::map<types::global_dof_index,TrilinosScalar>::const_iterator
287  dof = boundary_values.begin();
288  dof != boundary_values.end();
289  ++dof)
290  if ((dof->first >= local_range.first) &&
291  (dof->first < local_range.second))
292  constrained_rows.push_back (dof->first);
293 
294  // then eliminate these rows and
295  // set their diagonal entry to
296  // what we have determined
297  // above. if the value already is
298  // nonzero, it will be preserved,
299  // in accordance with the basic
300  // matrix classes in deal.II.
301  matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
302 
303  std::vector<types::global_dof_index> indices;
304  std::vector<TrilinosScalar> solution_values;
305  for (std::map<types::global_dof_index,TrilinosScalar>::const_iterator
306  dof = boundary_values.begin();
307  dof != boundary_values.end();
308  ++dof)
309  if ((dof->first >= local_range.first) &&
310  (dof->first < local_range.second))
311  {
312  indices.push_back (dof->first);
313  solution_values.push_back (dof->second);
314  }
315  solution.set (indices, solution_values);
316 
317  // now also set appropriate
318  // values for the rhs
319  for (unsigned int i=0; i<solution_values.size(); ++i)
320  solution_values[i] *= matrix.diag_element(indices[i]);
321 
322  right_hand_side.set (indices, solution_values);
323  }
324  else
325  {
326  // clear_rows() is a collective operation so we still have to call
327  // it:
328  std::vector<types::global_dof_index> constrained_rows;
329  matrix.clear_rows (constrained_rows, 1.);
330  }
331 
332  // clean up
333  matrix.compress (VectorOperation::insert);
334  solution.compress (VectorOperation::insert);
335  right_hand_side.compress (VectorOperation::insert);
336  }
337 
338 
339 
340  template <typename TrilinosMatrix, typename TrilinosBlockVector>
341  void
342  apply_block_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
343  TrilinosMatrix &matrix,
344  TrilinosBlockVector &solution,
345  TrilinosBlockVector &right_hand_side,
346  const bool eliminate_columns)
347  {
348  Assert (eliminate_columns == false, ExcNotImplemented());
349 
350  Assert (matrix.n() == right_hand_side.size(),
351  ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
352  Assert (matrix.n() == solution.size(),
353  ExcDimensionMismatch(matrix.n(), solution.size()));
354  Assert (matrix.n_block_rows() == matrix.n_block_cols(),
355  ExcNotQuadratic());
356 
357  const unsigned int n_blocks = matrix.n_block_rows();
358 
359  // We need to find the subdivision
360  // into blocks for the boundary values.
361  // To this end, generate a vector of
362  // maps with the respective indices.
363  std::vector<std::map<types::global_dof_index,TrilinosScalar> > block_boundary_values(n_blocks);
364  {
365  int block=0;
366  types::global_dof_index offset = 0;
367  for (std::map<types::global_dof_index,TrilinosScalar>::const_iterator
368  dof = boundary_values.begin();
369  dof != boundary_values.end();
370  ++dof)
371  {
372  if (dof->first >= matrix.block(block,0).m() + offset)
373  {
374  offset += matrix.block(block,0).m();
375  block++;
376  }
377  const types::global_dof_index index = dof->first - offset;
378  block_boundary_values[block].insert(
379  std::pair<types::global_dof_index, TrilinosScalar> (index,dof->second));
380  }
381  }
382 
383  // Now call the non-block variants on
384  // the diagonal subblocks and the
385  // solution/rhs.
386  for (unsigned int block=0; block<n_blocks; ++block)
387  TrilinosWrappers::apply_boundary_values(block_boundary_values[block],
388  matrix.block(block,block),
389  solution.block(block),
390  right_hand_side.block(block),
391  eliminate_columns);
392 
393  // Finally, we need to do something
394  // about the off-diagonal matrices. This
395  // is luckily not difficult. Just clear
396  // the whole row.
397  for (unsigned int block_m=0; block_m<n_blocks; ++block_m)
398  {
399  const std::pair<types::global_dof_index, types::global_dof_index> local_range
400  = matrix.block(block_m,0).local_range();
401 
402  std::vector<types::global_dof_index> constrained_rows;
403  for (std::map<types::global_dof_index,TrilinosScalar>::const_iterator
404  dof = block_boundary_values[block_m].begin();
405  dof != block_boundary_values[block_m].end();
406  ++dof)
407  if ((dof->first >= local_range.first) &&
408  (dof->first < local_range.second))
409  constrained_rows.push_back (dof->first);
410 
411  for (unsigned int block_n=0; block_n<n_blocks; ++block_n)
412  if (block_m != block_n)
413  matrix.block(block_m,block_n).clear_rows(constrained_rows);
414  }
415  }
416  }
417  }
418 
419 
420 
421  void
422  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
425  TrilinosWrappers::MPI::Vector &right_hand_side,
426  const bool eliminate_columns)
427  {
428  // simply redirect to the generic function
429  // used for both trilinos matrix types
430  internal::TrilinosWrappers::apply_boundary_values (boundary_values, matrix, solution,
431  right_hand_side, eliminate_columns);
432  }
433 
434 
435 
436  void
437  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
440  TrilinosWrappers::MPI::BlockVector &right_hand_side,
441  const bool eliminate_columns)
442  {
443  internal::TrilinosWrappers::apply_block_boundary_values (boundary_values, matrix,
444  solution, right_hand_side,
445  eliminate_columns);
446  }
447 
448 #endif
449 
450 } // namespace MatrixTools
451 
452 DEAL_II_NAMESPACE_CLOSE
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
Contents is actually a matrix.
std::pair< size_type, size_type > local_range() const
std::size_t size() const
void compress(const VectorOperation::values operation)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcNotImplemented()
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcInternalError()
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)