38#ifdef DEAL_II_WITH_PETSC
45#ifdef DEAL_II_WITH_TRILINOS
60#ifdef DEAL_II_WITH_PETSC
63 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
67 const bool eliminate_columns)
77 if (boundary_values.size() > 0)
79 const std::pair<types::global_dof_index, types::global_dof_index>
80 local_range = matrix.local_range();
89 PetscScalar average_nonzero_diagonal_entry = 1;
91 i < local_range.second;
93 if (matrix.diag_element(i) != PetscScalar())
95 average_nonzero_diagonal_entry =
std::abs(matrix.diag_element(i));
101 std::vector<types::global_dof_index> constrained_rows;
102 for (
const auto &boundary_value : boundary_values)
103 if ((boundary_value.first >= local_range.first) &&
104 (boundary_value.first < local_range.second))
105 constrained_rows.push_back(boundary_value.first);
118 if (eliminate_columns)
119 matrix.clear_rows_columns(constrained_rows,
120 average_nonzero_diagonal_entry);
122 matrix.clear_rows(constrained_rows, average_nonzero_diagonal_entry);
124 std::vector<types::global_dof_index> indices;
125 std::vector<PetscScalar> solution_values;
126 for (
const auto &boundary_value : boundary_values)
127 if ((boundary_value.first >= local_range.first) &&
128 (boundary_value.first < local_range.second))
130 indices.push_back(boundary_value.first);
131 solution_values.push_back(boundary_value.second);
133 solution.
set(indices, solution_values);
136 for (
auto &solution_value : solution_values)
137 solution_value *= average_nonzero_diagonal_entry;
139 right_hand_side.
set(indices, solution_values);
145 std::vector<types::global_dof_index> constrained_rows;
146 if (eliminate_columns)
147 matrix.clear_rows_columns(constrained_rows, 1.);
149 matrix.clear_rows(constrained_rows, 1.);
160 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
164 const bool eliminate_columns)
173 const unsigned int n_blocks = matrix.n_block_rows();
179 std::vector<std::map<::types::global_dof_index, PetscScalar>>
180 block_boundary_values(n_blocks);
184 for (
const auto &boundary_value : boundary_values)
186 if (boundary_value.first >= matrix.block(block, 0).m() + offset)
188 offset += matrix.block(block, 0).m();
192 block_boundary_values[block].insert(
193 std::pair<types::global_dof_index, PetscScalar>(
194 index, boundary_value.second));
201 for (
unsigned int block = 0; block < n_blocks; ++block)
203 matrix.block(block, block),
204 solution.
block(block),
205 right_hand_side.
block(block),
212 for (
unsigned int block_m = 0; block_m < n_blocks; ++block_m)
214 const std::pair<types::global_dof_index, types::global_dof_index>
215 local_range = matrix.block(block_m, 0).local_range();
217 std::vector<types::global_dof_index> constrained_rows;
218 for (std::map<types::global_dof_index, PetscScalar>::const_iterator
219 dof = block_boundary_values[block_m].begin();
220 dof != block_boundary_values[block_m].end();
222 if ((dof->first >= local_range.first) &&
223 (dof->first < local_range.second))
224 constrained_rows.push_back(dof->first);
226 for (
unsigned int block_n = 0; block_n < n_blocks; ++block_n)
227 if (block_m != block_n)
228 matrix.block(block_m, block_n).clear_rows(constrained_rows);
236#ifdef DEAL_II_WITH_TRILINOS
242 template <
typename TrilinosMatrix,
typename TrilinosVector>
246 TrilinosMatrix & matrix,
247 TrilinosVector & solution,
248 TrilinosVector & right_hand_side,
249 const bool eliminate_columns)
252 (void)eliminate_columns;
254 Assert(matrix.n() == right_hand_side.size(),
256 Assert(matrix.n() == solution.size(),
262 if (boundary_values.size() > 0)
264 const std::pair<types::global_dof_index, types::global_dof_index>
265 local_range = matrix.local_range();
266 Assert(local_range == right_hand_side.local_range(),
276 i < local_range.second;
278 if (matrix.diag_element(i) != 0)
280 average_nonzero_diagonal_entry =
281 std::fabs(matrix.diag_element(i));
287 std::vector<types::global_dof_index> constrained_rows;
288 for (
const auto &boundary_value : boundary_values)
289 if ((boundary_value.first >= local_range.first) &&
290 (boundary_value.first < local_range.second))
291 constrained_rows.push_back(boundary_value.first);
300 matrix.clear_rows(constrained_rows, average_nonzero_diagonal_entry);
302 std::vector<types::global_dof_index> indices;
303 std::vector<TrilinosScalar> solution_values;
304 for (
const auto &boundary_value : boundary_values)
305 if ((boundary_value.first >= local_range.first) &&
306 (boundary_value.first < local_range.second))
308 indices.push_back(boundary_value.first);
309 solution_values.push_back(boundary_value.second);
311 solution.set(indices, solution_values);
315 for (
unsigned int i = 0; i < solution_values.size(); ++i)
316 solution_values[i] *= matrix.diag_element(indices[i]);
318 right_hand_side.set(indices, solution_values);
324 std::vector<types::global_dof_index> constrained_rows;
325 matrix.clear_rows(constrained_rows, 1.);
336 template <
typename TrilinosMatrix,
typename TrilinosBlockVector>
339 const std::map<types::global_dof_index, TrilinosScalar>
341 TrilinosMatrix & matrix,
342 TrilinosBlockVector &solution,
343 TrilinosBlockVector &right_hand_side,
344 const bool eliminate_columns)
348 Assert(matrix.n() == right_hand_side.size(),
350 Assert(matrix.n() == solution.size(),
352 Assert(matrix.n_block_rows() == matrix.n_block_cols(),
355 const unsigned int n_blocks = matrix.n_block_rows();
361 std::vector<std::map<types::global_dof_index, TrilinosScalar>>
362 block_boundary_values(n_blocks);
366 for (
const auto &boundary_value : boundary_values)
368 if (boundary_value.first >= matrix.block(block, 0).m() + offset)
370 offset += matrix.block(block, 0).m();
374 boundary_value.first - offset;
375 block_boundary_values[block].insert(
376 std::pair<types::global_dof_index, TrilinosScalar>(
377 index, boundary_value.second));
384 for (
unsigned int block = 0; block < n_blocks; ++block)
386 matrix.block(block, block),
387 solution.block(block),
388 right_hand_side.block(block),
395 for (
unsigned int block_m = 0; block_m < n_blocks; ++block_m)
397 const std::pair<types::global_dof_index, types::global_dof_index>
398 local_range = matrix.block(block_m, 0).local_range();
400 std::vector<types::global_dof_index> constrained_rows;
403 block_boundary_values[block_m].begin();
404 dof != block_boundary_values[block_m].end();
406 if ((dof->first >= local_range.first) &&
407 (dof->first < local_range.second))
408 constrained_rows.push_back(dof->first);
410 for (
unsigned int block_n = 0; block_n < n_blocks; ++block_n)
411 if (block_m != block_n)
412 matrix.block(block_m, block_n).clear_rows(constrained_rows);
422 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
426 const bool eliminate_columns)
431 boundary_values, matrix, solution, right_hand_side, eliminate_columns);
438 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
442 const bool eliminate_columns)
445 boundary_values, matrix, solution, right_hand_side, eliminate_columns);
BlockType & block(const unsigned int i)
std::pair< size_type, size_type > local_range() const
void compress(const VectorOperation::values operation)
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)