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_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> 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> 55 DEAL_II_NAMESPACE_OPEN
60 #ifdef DEAL_II_WITH_PETSC 63 (
const std::map<types::global_dof_index,PetscScalar> &boundary_values,
67 const bool eliminate_columns)
69 (void)eliminate_columns;
80 if (boundary_values.size() > 0)
82 const std::pair<types::global_dof_index, types::global_dof_index> local_range
83 = matrix.local_range();
93 PetscScalar average_nonzero_diagonal_entry = 1;
95 if (matrix.diag_element(i) != PetscScalar ())
97 average_nonzero_diagonal_entry = std::abs(matrix.diag_element(i));
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();
108 if ((dof->first >= local_range.first) &&
109 (dof->first < local_range.second))
110 constrained_rows.push_back (dof->first);
123 matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
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();
131 if ((dof->first >= local_range.first) &&
132 (dof->first < local_range.second))
134 indices.push_back (dof->first);
135 solution_values.push_back (dof->second);
137 solution.
set (indices, solution_values);
141 for (
unsigned int i=0; i<solution_values.size(); ++i)
142 solution_values[i] *= average_nonzero_diagonal_entry;
144 right_hand_side.
set (indices, solution_values);
150 std::vector<types::global_dof_index> constrained_rows;
151 matrix.clear_rows (constrained_rows, 1.);
165 const bool eliminate_columns)
171 Assert (matrix.n_block_rows() == matrix.n_block_cols(),
174 const unsigned int n_blocks = matrix.n_block_rows();
180 std::vector<std::map<::types::global_dof_index,PetscScalar> > block_boundary_values(n_blocks);
184 for (std::map<types::global_dof_index,PetscScalar>::const_iterator
185 dof = boundary_values.begin();
186 dof != boundary_values.end();
189 if (dof->first >= matrix.block(block,0).m() + offset)
191 offset += matrix.block(block,0).m();
195 block_boundary_values[block].insert(std::pair<types::global_dof_index, PetscScalar> (index,dof->second));
202 for (
unsigned int block=0; block<n_blocks; ++block)
204 matrix.block(block,block),
205 solution.
block(block),
206 right_hand_side.
block(block),
213 for (
unsigned int block_m=0; block_m<n_blocks; ++block_m)
215 const std::pair<types::global_dof_index, types::global_dof_index> local_range
216 = matrix.block(block_m,0).local_range();
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();
223 if ((dof->first >= local_range.first) &&
224 (dof->first < local_range.second))
225 constrained_rows.push_back (dof->first);
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);
237 #ifdef DEAL_II_WITH_TRILINOS 243 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> local_range
265 = matrix.local_range();
266 Assert (local_range == right_hand_side.local_range(),
268 Assert (local_range == solution.local_range(),
275 TrilinosScalar average_nonzero_diagonal_entry = 1;
277 if (matrix.diag_element(i) != 0)
279 average_nonzero_diagonal_entry = std::fabs(matrix.diag_element(i));
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();
290 if ((dof->first >= local_range.first) &&
291 (dof->first < local_range.second))
292 constrained_rows.push_back (dof->first);
301 matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry);
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();
309 if ((dof->first >= local_range.first) &&
310 (dof->first < local_range.second))
312 indices.push_back (dof->first);
313 solution_values.push_back (dof->second);
315 solution.set (indices, solution_values);
319 for (
unsigned int i=0; i<solution_values.size(); ++i)
320 solution_values[i] *=
matrix.diag_element(indices[i]);
322 right_hand_side.set (indices, solution_values);
328 std::vector<types::global_dof_index> constrained_rows;
329 matrix.clear_rows (constrained_rows, 1.);
340 template <
typename TrilinosMatrix,
typename TrilinosBlockVector>
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)
357 const unsigned int n_blocks =
matrix.n_block_rows();
363 std::vector<std::map<types::global_dof_index,TrilinosScalar> > block_boundary_values(n_blocks);
367 for (std::map<types::global_dof_index,TrilinosScalar>::const_iterator
368 dof = boundary_values.begin();
369 dof != boundary_values.end();
372 if (dof->first >=
matrix.block(block,0).m() + offset)
374 offset +=
matrix.block(block,0).m();
378 block_boundary_values[block].insert(
379 std::pair<types::global_dof_index, TrilinosScalar> (index,dof->second));
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),
397 for (
unsigned int block_m=0; block_m<n_blocks; ++block_m)
399 const std::pair<types::global_dof_index, types::global_dof_index> local_range
400 =
matrix.block(block_m,0).local_range();
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();
407 if ((dof->first >= local_range.first) &&
408 (dof->first < local_range.second))
409 constrained_rows.push_back (dof->first);
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);
426 const bool eliminate_columns)
430 internal::TrilinosWrappers::apply_boundary_values (boundary_values, matrix, solution,
431 right_hand_side, eliminate_columns);
441 const bool eliminate_columns)
443 internal::TrilinosWrappers::apply_block_boundary_values (boundary_values, matrix,
444 solution, right_hand_side,
452 DEAL_II_NAMESPACE_CLOSE
Contents is actually a matrix.
std::pair< size_type, size_type > local_range() const
void compress(const VectorOperation::values operation)
unsigned int global_dof_index
#define Assert(cond, exc)
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)