42 #ifdef DEAL_II_WITH_PETSC
49 #ifdef DEAL_II_WITH_TRILINOS
64 #ifdef DEAL_II_WITH_PETSC
67 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
71 const bool eliminate_columns)
73 (void)eliminate_columns;
84 if (boundary_values.size() > 0)
86 const std::pair<types::global_dof_index, types::global_dof_index>
87 local_range =
matrix.local_range();
96 PetscScalar average_nonzero_diagonal_entry = 1;
98 i < local_range.second;
100 if (
matrix.diag_element(i) != PetscScalar())
102 average_nonzero_diagonal_entry = std::abs(
matrix.diag_element(i));
108 std::vector<types::global_dof_index> constrained_rows;
109 for (
const auto &boundary_value : boundary_values)
110 if ((boundary_value.first >= local_range.first) &&
111 (boundary_value.first < local_range.second))
112 constrained_rows.push_back(boundary_value.first);
125 matrix.clear_rows(constrained_rows, average_nonzero_diagonal_entry);
127 std::vector<types::global_dof_index> indices;
128 std::vector<PetscScalar> solution_values;
129 for (
const auto &boundary_value : boundary_values)
130 if ((boundary_value.first >= local_range.first) &&
131 (boundary_value.first < local_range.second))
133 indices.push_back(boundary_value.first);
134 solution_values.push_back(boundary_value.second);
136 solution.
set(indices, solution_values);
139 for (
auto &solution_value : solution_values)
140 solution_value *= average_nonzero_diagonal_entry;
142 right_hand_side.
set(indices, solution_values);
148 std::vector<types::global_dof_index> constrained_rows;
149 matrix.clear_rows(constrained_rows, 1.);
160 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
164 const bool eliminate_columns)
178 std::vector<std::map<::types::global_dof_index, PetscScalar>>
183 for (
const auto &boundary_value : boundary_values)
185 if (boundary_value.first >=
matrix.block(block, 0).m() + offset)
187 offset +=
matrix.block(block, 0).m();
191 block_boundary_values[block].insert(
192 std::pair<types::global_dof_index, PetscScalar>(
193 index, boundary_value.second));
200 for (
unsigned int block = 0; block <
n_blocks; ++block)
202 matrix.block(block, block),
203 solution.
block(block),
204 right_hand_side.
block(block),
211 for (
unsigned int block_m = 0; block_m <
n_blocks; ++block_m)
213 const std::pair<types::global_dof_index, types::global_dof_index>
214 local_range =
matrix.block(block_m, 0).local_range();
216 std::vector<types::global_dof_index> constrained_rows;
217 for (std::map<types::global_dof_index, PetscScalar>::const_iterator
218 dof = block_boundary_values[block_m].
begin();
219 dof != block_boundary_values[block_m].end();
221 if ((dof->first >= local_range.first) &&
222 (dof->first < local_range.second))
223 constrained_rows.push_back(dof->first);
225 for (
unsigned int block_n = 0; block_n <
n_blocks; ++block_n)
226 if (block_m != block_n)
227 matrix.block(block_m, block_n).clear_rows(constrained_rows);
235 #ifdef DEAL_II_WITH_TRILINOS
241 template <
typename TrilinosMatrix,
typename TrilinosVector>
246 TrilinosVector & solution,
247 TrilinosVector & right_hand_side,
248 const bool eliminate_columns)
251 (void)eliminate_columns;
261 if (boundary_values.size() > 0)
263 const std::pair<types::global_dof_index, types::global_dof_index>
264 local_range =
matrix.local_range();
265 Assert(local_range == right_hand_side.local_range(),
275 i < local_range.second;
277 if (
matrix.diag_element(i) != 0)
279 average_nonzero_diagonal_entry =
286 std::vector<types::global_dof_index> constrained_rows;
287 for (
const auto &boundary_value : boundary_values)
288 if ((boundary_value.first >= local_range.first) &&
289 (boundary_value.first < local_range.second))
290 constrained_rows.push_back(boundary_value.first);
299 matrix.clear_rows(constrained_rows, average_nonzero_diagonal_entry);
301 std::vector<types::global_dof_index> indices;
302 std::vector<TrilinosScalar> solution_values;
303 for (
const auto &boundary_value : boundary_values)
304 if ((boundary_value.first >= local_range.first) &&
305 (boundary_value.first < local_range.second))
307 indices.push_back(boundary_value.first);
308 solution_values.push_back(boundary_value.second);
310 solution.set(indices, solution_values);
314 for (
unsigned int i = 0; i < solution_values.size(); ++i)
315 solution_values[i] *=
matrix.diag_element(indices[i]);
317 right_hand_side.set(indices, solution_values);
323 std::vector<types::global_dof_index> constrained_rows;
324 matrix.clear_rows(constrained_rows, 1.);
335 template <
typename TrilinosMatrix,
typename TrilinosBlockVector>
338 const std::map<types::global_dof_index, TrilinosScalar>
341 TrilinosBlockVector &solution,
342 TrilinosBlockVector &right_hand_side,
343 const bool eliminate_columns)
360 std::vector<std::map<types::global_dof_index, TrilinosScalar>>
365 for (
const auto &boundary_value : boundary_values)
367 if (boundary_value.first >=
matrix.block(block, 0).m() + offset)
369 offset +=
matrix.block(block, 0).m();
373 boundary_value.first - offset;
374 block_boundary_values[block].insert(
375 std::pair<types::global_dof_index, TrilinosScalar>(
376 index, boundary_value.second));
383 for (
unsigned int block = 0; block <
n_blocks; ++block)
385 matrix.block(block, block),
386 solution.block(block),
387 right_hand_side.block(block),
394 for (
unsigned int block_m = 0; block_m <
n_blocks; ++block_m)
396 const std::pair<types::global_dof_index, types::global_dof_index>
397 local_range =
matrix.block(block_m, 0).local_range();
399 std::vector<types::global_dof_index> constrained_rows;
402 block_boundary_values[block_m].
begin();
403 dof != block_boundary_values[block_m].end();
405 if ((dof->first >= local_range.first) &&
406 (dof->first < local_range.second))
407 constrained_rows.push_back(dof->first);
409 for (
unsigned int block_n = 0; block_n <
n_blocks; ++block_n)
410 if (block_m != block_n)
411 matrix.block(block_m, block_n).clear_rows(constrained_rows);
421 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
425 const bool eliminate_columns)
430 boundary_values,
matrix, solution, right_hand_side, eliminate_columns);
437 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
441 const bool eliminate_columns)
444 boundary_values,
matrix, solution, right_hand_side, eliminate_columns);