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> 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_vector.h> 42 # include <deal.II/lac/petsc_parallel_block_sparse_matrix.h> 45 #ifdef DEAL_II_WITH_TRILINOS 46 # include <deal.II/lac/trilinos_sparse_matrix.h> 47 # include <deal.II/lac/trilinos_vector.h> 48 # include <deal.II/lac/trilinos_block_sparse_matrix.h> 49 # include <deal.II/lac/trilinos_block_vector.h> 57 DEAL_II_NAMESPACE_OPEN
62 #ifdef DEAL_II_WITH_PETSC 68 template <
typename PETScMatrix,
typename PETScVector>
72 PETScVector &solution,
73 PETScVector &right_hand_side,
74 const bool eliminate_columns)
76 (void)eliminate_columns;
77 Assert (eliminate_columns ==
false, ExcNotImplemented());
79 Assert (matrix.n() == right_hand_side.size(),
80 ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
81 Assert (matrix.n() == solution.size(),
82 ExcDimensionMismatch(matrix.n(), solution.size()));
87 if (boundary_values.size() > 0)
89 const std::pair<types::global_dof_index, types::global_dof_index> local_range
90 = matrix.local_range();
91 Assert (local_range == right_hand_side.local_range(),
93 Assert (local_range == solution.local_range(),
100 PetscScalar average_nonzero_diagonal_entry = 1;
102 if (matrix.diag_element(i) != PetscScalar ())
104 average_nonzero_diagonal_entry = std::abs(matrix.diag_element(i));
110 std::vector<types::global_dof_index> constrained_rows;
111 for (std::map<types::global_dof_index,double>::const_iterator
112 dof = boundary_values.begin();
113 dof != boundary_values.end();
115 if ((dof->first >= local_range.first) &&
116 (dof->first < local_range.second))
117 constrained_rows.push_back (dof->first);
130 matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
132 std::vector<types::global_dof_index> indices;
133 std::vector<PetscScalar> solution_values;
134 for (std::map<types::global_dof_index,double>::const_iterator
135 dof = boundary_values.begin();
136 dof != boundary_values.end();
138 if ((dof->first >= local_range.first) &&
139 (dof->first < local_range.second))
141 indices.push_back (dof->first);
142 solution_values.push_back (dof->second);
144 solution.set (indices, solution_values);
148 for (
unsigned int i=0; i<solution_values.size(); ++i)
149 solution_values[i] *= average_nonzero_diagonal_entry;
151 right_hand_side.set (indices, solution_values);
157 std::vector<types::global_dof_index> constrained_rows;
158 matrix.clear_rows (constrained_rows, 1.);
162 solution.compress (VectorOperation::insert);
163 right_hand_side.compress (VectorOperation::insert);
177 const bool eliminate_columns)
181 internal::PETScWrappers::apply_boundary_values (boundary_values, matrix, solution,
182 right_hand_side, eliminate_columns);
192 const bool eliminate_columns)
196 internal::PETScWrappers::apply_boundary_values (boundary_values, matrix, solution,
197 right_hand_side, eliminate_columns);
206 const bool eliminate_columns)
209 ExcDimensionMismatch(matrix.
n(), right_hand_side.
size()));
211 ExcDimensionMismatch(matrix.
n(), solution.
size()));
221 std::vector<std::map<::types::global_dof_index,double> > block_boundary_values(n_blocks);
225 for (std::map<types::global_dof_index,double>::const_iterator
226 dof = boundary_values.begin();
227 dof != boundary_values.end();
230 if (dof->first >= matrix.
block(block,0).
m() + offset)
232 offset += matrix.
block(block,0).
m();
236 block_boundary_values[block].insert(std::pair<types::global_dof_index, double> (index,dof->second));
243 for (
unsigned int block=0; block<n_blocks; ++block)
244 internal::PETScWrappers::apply_boundary_values(block_boundary_values[block],
245 matrix.
block(block,block),
246 solution.
block(block),
247 right_hand_side.
block(block),
254 for (
unsigned int block_m=0; block_m<n_blocks; ++block_m)
256 const std::pair<types::global_dof_index, types::global_dof_index> local_range
257 = matrix.
block(block_m,0).local_range();
259 std::vector<types::global_dof_index> constrained_rows;
260 for (std::map<types::global_dof_index,double>::const_iterator
261 dof = block_boundary_values[block_m].begin();
262 dof != block_boundary_values[block_m].end();
264 if ((dof->first >= local_range.first) &&
265 (dof->first < local_range.second))
266 constrained_rows.push_back (dof->first);
268 for (
unsigned int block_n=0; block_n<n_blocks; ++block_n)
269 if (block_m != block_n)
270 matrix.
block(block_m,block_n).clear_rows(constrained_rows);
278 #ifdef DEAL_II_WITH_TRILINOS 284 template <
typename TrilinosMatrix,
typename TrilinosVector>
287 TrilinosMatrix &matrix,
288 TrilinosVector &solution,
289 TrilinosVector &right_hand_side,
290 const bool eliminate_columns)
292 Assert (eliminate_columns ==
false, ExcNotImplemented());
293 (void)eliminate_columns;
295 Assert (matrix.n() == right_hand_side.size(),
296 ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
297 Assert (matrix.n() == solution.size(),
298 ExcDimensionMismatch(matrix.m(), solution.size()));
303 if (boundary_values.size() > 0)
305 const std::pair<types::global_dof_index, types::global_dof_index> local_range
306 = matrix.local_range();
307 Assert (local_range == right_hand_side.local_range(),
309 Assert (local_range == solution.local_range(),
316 TrilinosScalar average_nonzero_diagonal_entry = 1;
318 if (matrix.diag_element(i) != 0)
320 average_nonzero_diagonal_entry = std::fabs(matrix.diag_element(i));
326 std::vector<types::global_dof_index> constrained_rows;
327 for (std::map<types::global_dof_index,double>::const_iterator
328 dof = boundary_values.begin();
329 dof != boundary_values.end();
331 if ((dof->first >= local_range.first) &&
332 (dof->first < local_range.second))
333 constrained_rows.push_back (dof->first);
342 matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
344 std::vector<types::global_dof_index> indices;
345 std::vector<TrilinosScalar> solution_values;
346 for (std::map<types::global_dof_index,double>::const_iterator
347 dof = boundary_values.begin();
348 dof != boundary_values.end();
350 if ((dof->first >= local_range.first) &&
351 (dof->first < local_range.second))
353 indices.push_back (dof->first);
354 solution_values.push_back (dof->second);
356 solution.set (indices, solution_values);
360 for (
unsigned int i=0; i<solution_values.size(); ++i)
361 solution_values[i] *= matrix.diag_element(indices[i]);
363 right_hand_side.set (indices, solution_values);
369 std::vector<types::global_dof_index> constrained_rows;
370 matrix.clear_rows (constrained_rows, 1.);
374 matrix.compress (VectorOperation::insert);
375 solution.compress (VectorOperation::insert);
376 right_hand_side.compress (VectorOperation::insert);
381 template <
typename TrilinosMatrix,
typename TrilinosBlockVector>
383 apply_block_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
384 TrilinosMatrix &matrix,
385 TrilinosBlockVector &solution,
386 TrilinosBlockVector &right_hand_side,
387 const bool eliminate_columns)
389 Assert (eliminate_columns ==
false, ExcNotImplemented());
391 Assert (matrix.n() == right_hand_side.size(),
392 ExcDimensionMismatch(matrix.n(), right_hand_side.size()));
393 Assert (matrix.n() == solution.size(),
394 ExcDimensionMismatch(matrix.n(), solution.size()));
395 Assert (matrix.n_block_rows() == matrix.n_block_cols(),
398 const unsigned int n_blocks = matrix.n_block_rows();
404 std::vector<std::map<types::global_dof_index,double> > block_boundary_values(n_blocks);
408 for (std::map<types::global_dof_index,double>::const_iterator
409 dof = boundary_values.begin();
410 dof != boundary_values.end();
413 if (dof->first >= matrix.block(block,0).m() + offset)
415 offset += matrix.block(block,0).m();
419 block_boundary_values[block].insert(
420 std::pair<types::global_dof_index, double> (index,dof->second));
427 for (
unsigned int block=0; block<n_blocks; ++block)
428 TrilinosWrappers::apply_boundary_values(block_boundary_values[block],
429 matrix.block(block,block),
430 solution.block(block),
431 right_hand_side.block(block),
438 for (
unsigned int block_m=0; block_m<n_blocks; ++block_m)
440 const std::pair<types::global_dof_index, types::global_dof_index> local_range
441 = matrix.block(block_m,0).local_range();
443 std::vector<types::global_dof_index> constrained_rows;
444 for (std::map<types::global_dof_index,double>::const_iterator
445 dof = block_boundary_values[block_m].begin();
446 dof != block_boundary_values[block_m].end();
448 if ((dof->first >= local_range.first) &&
449 (dof->first < local_range.second))
450 constrained_rows.push_back (dof->first);
452 for (
unsigned int block_n=0; block_n<n_blocks; ++block_n)
453 if (block_m != block_n)
454 matrix.block(block_m,block_n).clear_rows(constrained_rows);
468 const bool eliminate_columns)
472 internal::TrilinosWrappers::apply_boundary_values (boundary_values, matrix, solution,
473 right_hand_side, eliminate_columns);
483 const bool eliminate_columns)
487 internal::TrilinosWrappers::apply_boundary_values (boundary_values, matrix, solution,
488 right_hand_side, eliminate_columns);
498 const bool eliminate_columns)
500 internal::TrilinosWrappers::apply_block_boundary_values (boundary_values, matrix,
501 solution, right_hand_side,
512 const bool eliminate_columns)
514 internal::TrilinosWrappers::apply_block_boundary_values (boundary_values, matrix,
515 solution, right_hand_side,
523 DEAL_II_NAMESPACE_CLOSE
unsigned int n_block_cols() const
unsigned int global_dof_index
#define Assert(cond, exc)
BlockType & block(const unsigned int row, const unsigned int column)
unsigned int n_block_rows() const
BlockType & block(const unsigned int i)