38#ifdef DEAL_II_WITH_PETSC
44#ifdef DEAL_II_WITH_TRILINOS
63 template <
typename Iterator>
65 column_less_than(
const typename Iterator::value_type p,
66 const unsigned int column)
68 return (p.column() < column);
74 template <
typename number>
77 const std::map<types::global_dof_index, number> &boundary_values,
81 const bool eliminate_columns)
87 Assert(matrix.n() == matrix.m(),
92 if (boundary_values.size() == 0)
104 number first_nonzero_diagonal_entry = 1;
105 for (
unsigned int i = 0; i < n_dofs; ++i)
106 if (matrix.diag_element(i) != number())
108 first_nonzero_diagonal_entry = matrix.diag_element(i);
113 typename std::map<types::global_dof_index, number>::const_iterator
114 dof = boundary_values.begin(),
115 endd = boundary_values.end();
116 for (; dof != endd; ++dof)
126 matrix.begin(dof_number);
127 p != matrix.end(dof_number);
129 if (p->column() != dof_number)
146 if (matrix.diag_element(dof_number) != number())
148 new_rhs = dof->second * matrix.diag_element(dof_number);
149 right_hand_side(dof_number) = new_rhs;
153 matrix.set(dof_number, dof_number, first_nonzero_diagonal_entry);
154 new_rhs = dof->second * first_nonzero_diagonal_entry;
155 right_hand_side(dof_number) = new_rhs;
166 if (eliminate_columns)
171 const number diagonal_entry = matrix.diag_element(dof_number);
182 matrix.begin(dof_number) + 1;
183 q != matrix.end(dof_number);
193 const unsigned int column) =
211 Assert((p != matrix.end(row)) && (p->column() == dof_number),
213 "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."));
220 right_hand_side(row) -=
221 static_cast<number
>(p->value()) / diagonal_entry * new_rhs;
229 solution(dof_number) = dof->second;
235 template <
typename number>
238 const std::map<types::global_dof_index, number> &boundary_values,
242 const bool eliminate_columns)
244 const unsigned int blocks = matrix.n_block_rows();
251 Assert(matrix.get_sparsity_pattern().get_row_indices() ==
252 matrix.get_sparsity_pattern().get_column_indices(),
254 Assert(matrix.get_sparsity_pattern().get_column_indices() ==
257 Assert(matrix.get_sparsity_pattern().get_row_indices() ==
263 if (boundary_values.size() == 0)
275 number first_nonzero_diagonal_entry = 0;
276 for (
unsigned int diag_block = 0; diag_block < blocks; ++diag_block)
278 for (
unsigned int i = 0; i < matrix.block(diag_block, diag_block).n();
280 if (matrix.block(diag_block, diag_block).diag_element(i) != number{})
282 first_nonzero_diagonal_entry =
283 matrix.block(diag_block, diag_block).diag_element(i);
289 if (first_nonzero_diagonal_entry != number{})
294 if (first_nonzero_diagonal_entry == number{})
295 first_nonzero_diagonal_entry = 1;
298 typename std::map<types::global_dof_index, number>::const_iterator
299 dof = boundary_values.begin(),
300 endd = boundary_values.end();
302 matrix.get_sparsity_pattern();
312 for (; dof != endd; ++dof)
321 const std::pair<unsigned int, types::global_dof_index> block_index =
331 for (
unsigned int block_col = 0; block_col < blocks; ++block_col)
333 (block_col == block_index.first ?
334 matrix.block(block_index.first, block_col)
335 .begin(block_index.second) +
337 matrix.block(block_index.first, block_col)
338 .begin(block_index.second));
339 p != matrix.block(block_index.first, block_col)
340 .end(block_index.second);
358 if (matrix.block(block_index.first, block_index.first)
359 .diag_element(block_index.second) != number{})
361 dof->second * matrix.block(block_index.first, block_index.first)
362 .diag_element(block_index.second);
365 matrix.block(block_index.first, block_index.first)
366 .diag_element(block_index.second) = first_nonzero_diagonal_entry;
367 new_rhs = dof->second * first_nonzero_diagonal_entry;
369 right_hand_side.
block(block_index.first)(block_index.second) = new_rhs;
381 if (eliminate_columns)
386 const number diagonal_entry =
387 matrix.block(block_index.first, block_index.first)
388 .diag_element(block_index.second);
411 for (
unsigned int block_row = 0; block_row < blocks; ++block_row)
416 sparsity_pattern.
block(block_row, block_index.first);
419 matrix.block(block_row, block_index.first);
421 matrix.block(block_index.first, block_row);
427 (block_index.first == block_row ?
428 transpose_matrix.
begin(block_index.second) + 1 :
429 transpose_matrix.
begin(block_index.second));
430 q != transpose_matrix.
end(block_index.second);
443 const unsigned int column) =
452 if (this_matrix.
begin(row)->column() ==
454 p = this_matrix.
begin(row);
457 this_matrix.
end(row),
463 this_matrix.
end(row),
477 Assert((p->column() == block_index.second) &&
478 (p != this_matrix.
end(row)),
482 right_hand_side.
block(block_row)(row) -=
483 number(p->value()) / diagonal_entry * new_rhs;
492 solution.
block(block_index.first)(block_index.second) = dof->second;
498 template <
typename number>
501 const std::map<types::global_dof_index, number> &boundary_values,
502 const std::vector<types::global_dof_index> & local_dof_indices,
505 const bool eliminate_columns)
507 Assert(local_dof_indices.size() == local_matrix.
m(),
509 Assert(local_dof_indices.size() == local_matrix.
n(),
511 Assert(local_dof_indices.size() == local_rhs.
size(),
516 if (boundary_values.size() == 0)
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)
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())
552 for (
unsigned int j = 0; j < n_local_dofs; ++j)
554 local_matrix(i, j) = 0;
561 if (local_matrix(i, i) == number{})
565 if (average_diagonal == number{})
567 unsigned int nonzero_diagonals = 0;
568 for (
unsigned int k = 0; k < n_local_dofs; ++k)
569 if (local_matrix(k, k) != number{})
571 average_diagonal +=
std::abs(local_matrix(k, k));
574 if (nonzero_diagonals != 0)
575 average_diagonal /= nonzero_diagonals;
577 average_diagonal = 0;
583 if (average_diagonal == number{})
584 average_diagonal = 1.;
586 local_matrix(i, i) = average_diagonal;
589 local_matrix(i, i) =
std::abs(local_matrix(i, i));
593 local_rhs(i) = local_matrix(i, i) * boundary_value->second;
597 if (eliminate_columns ==
true)
599 for (
unsigned int row = 0; row < n_local_dofs; ++row)
603 local_matrix(row, i) * boundary_value->second;
604 local_matrix(row, i) = 0;
615#include "matrix_tools.inst"
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
SparsityPatternType & block(const size_type row, const size_type column)
const BlockIndices & get_column_indices() const
BlockType & block(const unsigned int i)
const BlockIndices & get_block_indices() const
const Accessor< number, Constness > & value_type
const_iterator end() const
const_iterator begin() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcBlocksDontMatch()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMessage(std::string arg1)
Iterator lower_bound(Iterator first, Iterator last, const T &val)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)