16 #include <deal.II/lac/constraint_matrix.h> 17 #include <deal.II/lac/constraint_matrix.templates.h> 19 #include <deal.II/base/memory_consumption.h> 20 #include <deal.II/lac/dynamic_sparsity_pattern.h> 21 #include <deal.II/lac/block_vector.h> 22 #include <deal.II/lac/block_sparse_matrix.h> 23 #include <deal.II/lac/sparse_matrix_ez.h> 24 #include <deal.II/lac/chunk_sparse_matrix.h> 25 #include <deal.II/lac/block_sparse_matrix_ez.h> 26 #include <deal.II/lac/la_vector.h> 27 #include <deal.II/lac/la_parallel_vector.h> 28 #include <deal.II/lac/la_parallel_block_vector.h> 29 #include <deal.II/lac/petsc_sparse_matrix.h> 30 #include <deal.II/lac/petsc_parallel_vector.h> 31 #include <deal.II/lac/petsc_parallel_block_vector.h> 32 #include <deal.II/lac/petsc_parallel_sparse_matrix.h> 33 #include <deal.II/lac/petsc_parallel_block_sparse_matrix.h> 34 #include <deal.II/lac/trilinos_vector.h> 35 #include <deal.II/lac/trilinos_parallel_block_vector.h> 36 #include <deal.II/lac/trilinos_sparse_matrix.h> 37 #include <deal.II/lac/trilinos_block_sparse_matrix.h> 38 #include <deal.II/lac/matrix_block.h> 39 #include <deal.II/lac/diagonal_matrix.h> 45 #include <boost/serialization/utility.hpp> 47 DEAL_II_NAMESPACE_OPEN
70 return (p.second == 0);
86 return index == a.
index;
103 return boost::make_iterator_range(
lines.begin(),
lines.end());
109 const IndexSet &locally_active_dofs,
110 const MPI_Comm mpi_communicator,
111 const bool verbose)
const 133 std::map< unsigned int, std::vector<ConstraintLine> > to_send;
135 const unsigned int myid = ::Utilities::MPI::this_mpi_process(mpi_communicator);
136 const unsigned int nproc = ::Utilities::MPI::n_mpi_processes(mpi_communicator);
140 IndexSet non_owned = locally_active_dofs;
142 for (
unsigned int owner=0; owner<nproc; ++owner)
145 IndexSet indices_to_send = non_owned & locally_owned_dofs[owner];
146 for (
const auto &row_idx : indices_to_send)
148 to_send[owner].push_back(get_line(row_idx));
154 unsigned int inconsistent = 0;
157 for (
const auto &kv : received)
160 for (
auto &lineit : kv.second)
169 std::cout <<
"Proc " << myid
170 <<
" got line " << lineit.index
171 <<
" from " << kv.first
172 <<
" inhomogeneity " << lineit.inhomogeneity <<
" != " << reference.
inhomogeneity << std::endl;
174 else if (lineit.entries != reference.
entries)
178 std::cout <<
"Proc " << myid
179 <<
" got line " << lineit.index
180 <<
" from " << kv.first
188 if (verbose && total>0 && myid==0)
189 std::cout << total <<
" inconsistent lines discovered!" << std::endl;
198 for (std::set<size_type>::const_iterator
209 if (
lines[i] ==
true)
227 const std::vector<std::pair<size_type,double> > &col_val_pairs)
241 col_val_pair = col_val_pairs.begin();
242 col_val_pair!=col_val_pairs.end(); ++col_val_pair)
244 Assert (line != col_val_pair->first,
245 ExcMessage (
"Can't constrain a degree of freedom to itself"));
247 for (ConstraintLine::Entries::const_iterator
249 p != line_ptr->
entries.end(); ++p)
250 if (p->first == col_val_pair->first)
253 Assert (p->second == col_val_pair->second,
255 p->second, col_val_pair->second));
259 line_ptr->
entries.push_back (*col_val_pair);
273 ExcMessage (
"Filter needs to be larger than constraint matrix size."));
274 for (std::vector<ConstraintLine>::const_iterator line=constraints.
lines.begin();
275 line!=constraints.
lines.end(); ++line)
281 for (
size_type i=0; i<line->entries.size(); ++i)
282 if (filter.
is_element(line->entries[i].first))
284 line->entries[i].second);
301 std::vector<size_type> new_lines (
lines_cache.size(),
304 for (std::vector<ConstraintLine>::const_iterator line=
lines.begin();
317 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
318 line!=
lines.end(); ++line)
322 line->entries.erase (std::remove_if (line->entries.begin(),
325 line->entries.end());
338 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
339 line!=
lines.end(); ++line)
341 for (ConstraintLine::Entries::iterator it = line->entries.begin(); it!=line->entries.end(); ++it)
343 largest_idx=std::max(largest_idx, it->first);
363 bool chained_constraint_replaced =
false;
365 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
366 line!=
lines.end(); ++line)
379 while (entry < line->entries.size())
387 chained_constraint_replaced =
true;
390 const size_type dof_index = line->entries[entry].first;
391 const double weight = line->entries[entry].second;
393 Assert (dof_index != line->index,
394 ExcMessage (
"Cycle in constraints detected!"));
409 if (constrained_line->
entries.size() > 0)
413 ExcMessage (
"Cycle in constraints detected!"));
417 line->entries[entry] =
418 std::make_pair (constrained_line->
entries[0].first,
419 constrained_line->
entries[0].second *
423 line->entries.emplace_back (constrained_line->
entries[i].first,
424 constrained_line->
entries[i].second
432 Assert(n_replacements/2<largest_idx,
ExcMessage(
"Cycle in constraints detected!"));
433 if (n_replacements/2>=largest_idx)
444 line->entries.erase (line->entries.begin()+entry);
462 if (chained_constraint_replaced ==
false)
476 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
477 line!=
lines.end(); ++line)
479 std::sort (line->entries.begin(), line->entries.end());
486 for (
size_type i=1; i<line->entries.size(); ++i)
487 if (line->entries[i].first == line->entries[i-1].first)
490 if (duplicates > 0 || line->entries.size() < line->entries.capacity())
497 new_entries = line->entries;
502 new_entries.reserve (line->entries.size() - duplicates);
503 new_entries.push_back(line->entries[0]);
504 for (
size_type j=1; j<line->entries.size(); ++j)
505 if (line->entries[j].first == line->entries[j-1].first)
507 Assert (new_entries.back().first == line->entries[j].first,
509 new_entries.back().second += line->entries[j].second;
512 new_entries.push_back (line->entries[j]);
514 Assert (new_entries.size() == line->entries.size() - duplicates,
519 for (
size_type j=1; j<new_entries.size(); ++j)
521 Assert (new_entries[j].first != new_entries[j-1].first,
523 Assert (new_entries[j].first > new_entries[j-1].first,
529 line->entries.swap (new_entries);
543 for (
size_type i=0; i<line->entries.size(); ++i)
544 sum += line->entries[i].second;
545 if ((sum != 1.0) && (std::fabs (sum-1.) < 1.e-13))
547 for (
size_type i=0; i<line->entries.size(); ++i)
548 line->entries[i].second /= sum;
549 line->inhomogeneity /= sum;
557 for (std::vector<ConstraintLine>::const_iterator line=
lines.begin();
558 line!=
lines.end(); ++line)
559 for (ConstraintLine::Entries::const_iterator
560 entry=line->entries.begin();
561 entry!=line->entries.end(); ++entry)
568 Assert (is_circle ==
false,
581 const bool allow_different_local_lines)
583 (void) allow_different_local_lines;
584 Assert(allow_different_local_lines ||
586 ExcMessage(
"local_lines for this and the other objects are not the same " 587 "although allow_different_local_lines is false."));
590 const bool object_was_sorted =
sorted;
602 for (std::vector<ConstraintLine>::iterator line=
lines.begin();
603 line!=
lines.end(); ++line)
606 for (
size_type i=0; i<line->entries.size(); ++i)
618 tmp.push_back(line->entries[i]);
626 Assert (other_line !=
nullptr,
629 const double weight = line->entries[i].second;
631 for (ConstraintLine::Entries::const_iterator j=other_line->begin();
632 j!=other_line->end(); ++j)
633 tmp.emplace_back(j->first, j->second*weight);
641 line->entries.swap (tmp);
654 for (std::vector<ConstraintLine>::const_iterator line =
lines.begin();
655 line !=
lines.end(); ++line)
664 for (std::vector<ConstraintLine>::const_iterator line = other_constraints.
lines.begin();
665 line != other_constraints.
lines.end(); ++line)
671 lines.push_back(*line);
677 lines.push_back(*line);
684 switch (merge_conflict_behavior)
715 if (object_was_sorted ==
true)
734 for (std::vector<ConstraintLine>::iterator i =
lines.begin();
735 i !=
lines.end(); ++i)
738 for (ConstraintLine::Entries::iterator
739 j = i->entries.begin();
740 j != i->entries.end(); ++j)
759 std::vector<ConstraintLine> tmp;
764 std::vector<size_type> tmp;
820 ((entry != sparsity.
end(row)) &&
821 entry->is_valid_entry());
824 const size_type column = entry->column();
845 const size_type column = entry->column();
923 sparsity.
add (row, new_col);
926 if ((new_col < column) && (old_rowlength != new_rowlength))
928 old_rowlength = new_rowlength;
951 .entries.size(); ++q)
991 const std::pair<size_type,size_type>
992 block_index = index_mapping.global_to_local(row);
993 const size_type block_row = block_index.first;
1002 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1005 block_sparsity = sparsity.
block(block_row, block_col);
1008 entry = block_sparsity.
begin(block_index.second);
1009 (entry != block_sparsity.
end(block_index.second)) &&
1014 = index_mapping.local_to_global(block_col, entry->column());
1032 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1035 block_sparsity = sparsity.
block(block_row,block_col);
1038 entry = block_sparsity.
begin(block_index.second);
1039 (entry != block_sparsity.
end(block_index.second)) &&
1044 = index_mapping.local_to_global (block_col, entry->column());
1099 for (
size_type row=0; row<n_rows; ++row)
1102 const std::pair<size_type,size_type>
1103 block_index = index_mapping.global_to_local(row);
1104 const size_type block_row = block_index.first;
1105 const size_type local_row = block_index.second;
1120 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1123 block_sparsity = sparsity.
block(block_row, block_col);
1128 = index_mapping.local_to_global(block_col,
1137 .entries.size(); ++q)
1148 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1151 block_sparsity = sparsity.
block(block_row,block_col);
1156 = index_mapping.local_to_global (block_col,
1196 return ((p.
entries.size() == 1) &&
1197 (p.
entries[0].second == 1.0));
1211 return ((p.
entries.size() == 1) &&
1212 (p.
entries[0].first == index2) &&
1213 (p.
entries[0].second == 1.0));
1222 return ((p.
entries.size() == 1) &&
1223 (p.
entries[0].first == index1) &&
1224 (p.
entries[0].second == 1.0));
1236 for (std::vector<ConstraintLine>::const_iterator i=
lines.begin();
1237 i!=
lines.end(); ++i)
1240 return_value = std::max(return_value,
1241 static_cast<size_type>(i->entries.size()));
1243 return return_value;
1250 for (std::vector<ConstraintLine>::const_iterator i=
lines.begin();
1251 i!=
lines.end(); ++i)
1252 if (i->inhomogeneity != 0.)
1264 if (
lines[i].entries.size() > 0)
1267 out <<
" " <<
lines[i].index
1268 <<
" " <<
lines[i].entries[j].first
1269 <<
": " <<
lines[i].entries[j].second <<
"\n";
1272 if (
lines[i].inhomogeneity != 0)
1273 out <<
" " <<
lines[i].index
1274 <<
": " <<
lines[i].inhomogeneity <<
"\n";
1281 if (
lines[i].inhomogeneity != 0)
1282 out <<
" " <<
lines[i].index
1283 <<
" = " <<
lines[i].inhomogeneity
1286 out <<
" " <<
lines[i].index <<
" = 0\n";
1298 out <<
"digraph constraints {" 1303 if (
lines[i].entries.size() > 0)
1305 out <<
" " <<
lines[i].index <<
"->" <<
lines[i].entries[j].first
1307 <<
lines[i].entries[j].second
1310 out <<
" " <<
lines[i].index <<
"\n";
1312 out <<
"}" << std::endl;
1331 const unsigned int indices_size = indices.size();
1332 const std::vector<std::pair<types::global_dof_index,double> > *line_ptr;
1333 for (
unsigned int i=0; i<indices_size; ++i)
1338 if (line_ptr!=
nullptr)
1340 const unsigned int line_size = line_ptr->size();
1341 for (
unsigned int j=0; j<line_size; ++j)
1342 indices.push_back((*line_ptr)[j].first);
1347 std::sort(indices.begin(),indices.end());
1348 std::vector<types::global_dof_index>::iterator it;
1349 it = std::unique(indices.begin(),indices.end());
1350 indices.resize(it-indices.begin());
1365 #define VECTOR_FUNCTIONS(VectorType) \ 1366 template void ConstraintMatrix::condense<VectorType >(const VectorType &uncondensed,\ 1367 VectorType &condensed) const;\ 1368 template void ConstraintMatrix::condense<VectorType >(VectorType &vec) const;\ 1369 template void ConstraintMatrix:: \ 1370 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1371 const std::vector<ConstraintMatrix::size_type> &, \ 1373 const FullMatrix<VectorType::value_type> &) const;\ 1374 template void ConstraintMatrix:: \ 1375 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1376 const std::vector<ConstraintMatrix::size_type> &, \ 1377 const std::vector<ConstraintMatrix::size_type> &, \ 1379 const FullMatrix<VectorType::value_type> &, \ 1382 #define PARALLEL_VECTOR_FUNCTIONS(VectorType) \ 1383 template void ConstraintMatrix:: \ 1384 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1385 const std::vector<ConstraintMatrix::size_type> &, \ 1387 const FullMatrix<VectorType::value_type> &) const;\ 1388 template void ConstraintMatrix:: \ 1389 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1390 const std::vector<ConstraintMatrix::size_type> &, \ 1391 const std::vector<ConstraintMatrix::size_type> &, \ 1393 const FullMatrix<VectorType::value_type> &, \ 1396 #ifdef DEAL_II_WITH_PETSC 1401 #ifdef DEAL_II_WITH_TRILINOS 1406 #define MATRIX_VECTOR_FUNCTIONS(MatrixType, VectorType) \ 1407 template void ConstraintMatrix:: \ 1408 distribute_local_to_global<MatrixType,VectorType > (const FullMatrix<MatrixType::value_type> &, \ 1409 const Vector<VectorType::value_type> &, \ 1410 const std::vector<ConstraintMatrix::size_type> &, \ 1414 std::integral_constant<bool, false>) const 1415 #define MATRIX_FUNCTIONS(MatrixType,VectorScalar) \ 1416 template void ConstraintMatrix:: \ 1417 distribute_local_to_global<MatrixType,Vector<VectorScalar> > (const FullMatrix<MatrixType::value_type> &, \ 1418 const Vector<VectorScalar> &, \ 1419 const std::vector<ConstraintMatrix::size_type> &, \ 1421 Vector<VectorScalar> &, \ 1423 std::integral_constant<bool, false>) const 1424 #define BLOCK_MATRIX_VECTOR_FUNCTIONS(MatrixType, VectorType) \ 1425 template void ConstraintMatrix:: \ 1426 distribute_local_to_global<MatrixType,VectorType > (const FullMatrix<MatrixType::value_type> &, \ 1427 const Vector<VectorType::value_type> &, \ 1428 const std::vector<ConstraintMatrix::size_type> &, \ 1432 std::integral_constant<bool, true>) const 1433 #define BLOCK_MATRIX_FUNCTIONS(MatrixType) \ 1434 template void ConstraintMatrix:: \ 1435 distribute_local_to_global<MatrixType,Vector<MatrixType::value_type> > (const FullMatrix<MatrixType::value_type> &, \ 1436 const Vector<MatrixType::value_type> &, \ 1437 const std::vector<ConstraintMatrix::size_type> &, \ 1439 Vector<MatrixType::value_type> &, \ 1441 std::integral_constant<bool, true>) const 1446 MATRIX_FUNCTIONS(
FullMatrix<std::complex<double> >,std::complex<double>);
1452 MATRIX_FUNCTIONS(
SparseMatrix<std::complex<double> >,std::complex<double>);
1453 MATRIX_FUNCTIONS(
SparseMatrix<std::complex<float> >,std::complex<float>);
1469 #ifdef DEAL_II_WITH_PETSC 1478 #ifdef DEAL_II_WITH_TRILINOS 1486 #define SPARSITY_FUNCTIONS(SparsityPatternType) \ 1487 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1488 const std::vector<ConstraintMatrix::size_type> &, \ 1489 SparsityPatternType &, \ 1491 const Table<2,bool> &, \ 1492 std::integral_constant<bool, false>) const; \ 1493 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1494 const std::vector<ConstraintMatrix::size_type> &, \ 1495 const std::vector<ConstraintMatrix::size_type> &, \ 1496 SparsityPatternType &, \ 1498 const Table<2,bool> &) const 1499 #define BLOCK_SPARSITY_FUNCTIONS(SparsityPatternType) \ 1500 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1501 const std::vector<ConstraintMatrix::size_type> &, \ 1502 SparsityPatternType &, \ 1504 const Table<2,bool> &, \ 1505 std::integral_constant<bool, true>) const; \ 1506 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1507 const std::vector<ConstraintMatrix::size_type> &, \ 1508 const std::vector<ConstraintMatrix::size_type> &, \ 1509 SparsityPatternType &, \ 1511 const Table<2,bool> &) const 1518 #ifdef DEAL_II_WITH_TRILINOS 1524 #define ONLY_MATRIX_FUNCTIONS(MatrixType) \ 1525 template void ConstraintMatrix::distribute_local_to_global<MatrixType > ( \ 1526 const FullMatrix<MatrixType::value_type> &, \ 1527 const std::vector<ConstraintMatrix::size_type> &, \ 1528 const std::vector<ConstraintMatrix::size_type> &, \ 1529 MatrixType &) const; \ 1530 template void ConstraintMatrix::distribute_local_to_global<MatrixType > ( \ 1531 const FullMatrix<MatrixType::value_type> &, \ 1532 const std::vector<ConstraintMatrix::size_type> &, \ 1533 const ConstraintMatrix &, \ 1534 const std::vector<ConstraintMatrix::size_type> &, \ 1546 #ifdef DEAL_II_WITH_TRILINOS 1551 #ifdef DEAL_II_WITH_PETSC 1557 #include "constraint_matrix.inst" 1565 #define SCRATCH_INITIALIZER(MatrixScalar,VectorScalar,Name) \ 1566 ConstraintMatrixData<MatrixScalar,VectorScalar>::ScratchData scratch_data_initializer_##Name; \ 1567 template <> Threads::ThreadLocalStorage<ConstraintMatrixData<MatrixScalar,VectorScalar>::ScratchData> \ 1568 ConstraintMatrixData<MatrixScalar,VectorScalar>::scratch_data(scratch_data_initializer_##Name) 1570 SCRATCH_INITIALIZER(
double,
double,dd);
1571 SCRATCH_INITIALIZER(
float,
float,ff);
1572 SCRATCH_INITIALIZER(std::complex<double>,std::complex<double>,zz);
1573 SCRATCH_INITIALIZER(std::complex<float>,std::complex<float>,cc);
1574 SCRATCH_INITIALIZER(
double,std::complex<double>,dz);
1575 SCRATCH_INITIALIZER(
float,std::complex<float>,fc);
1576 #undef SCRATCH_INITIALIZER 1580 DEAL_II_NAMESPACE_CLOSE
const types::global_dof_index invalid_size_type
bool operator<(const ConstraintLine &) const
std::vector< size_type > lines_cache
static ::ExceptionBase & ExcDoFIsConstrainedFromBothObjects(size_type arg1)
static ::ExceptionBase & ExcEntryAlreadyExists(size_type arg1, size_type arg2, double arg3, double arg4)
std::vector< ConstraintLine >::const_iterator const_iterator
void merge(const ConstraintMatrix &other_constraints, const MergeConflictBehavior merge_conflict_behavior=no_conflicts_allowed, const bool allow_different_local_lines=false)
static bool check_zero_weight(const std::pair< size_type, double > &p)
static ::ExceptionBase & ExcIO()
std::size_t memory_consumption() const
size_type column_number(const size_type row, const size_type index) const
void add(const size_type i, const size_type j)
size_type n_block_cols() const
boost::iterator_range< const_iterator > LineRange
size_type n_block_rows() const
#define AssertIndexRange(index, range)
void add_entries(const size_type line, const std::vector< std::pair< size_type, double > > &col_val_pairs)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
#define AssertThrow(cond, exc)
types::global_dof_index size_type
std::vector< std::pair< size_type, double > > Entries
static ::ExceptionBase & ExcLineInexistant(size_type arg1)
bool is_compressed() const
SparsityPatternType & block(const size_type row, const size_type column)
bool is_valid_entry() const
void add_line(const size_type line)
bool operator==(const ConstraintLine &) const
static ::ExceptionBase & ExcMessage(std::string arg1)
void add_entry(const size_type line, const size_type column, const double value)
void add(const size_type i, const size_type j)
void add(const size_type i, const size_type j)
T sum(const T &t, const MPI_Comm &mpi_communicator)
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
void print(std::ostream &out) const
size_type index_within_set(const size_type global_index) const
void condense(SparsityPattern &sparsity) const
size_type n_constraints() const
const BlockIndices & get_row_indices() const
bool is_consistent_in_parallel(const std::vector< IndexSet > &locally_owned_dofs, const IndexSet &locally_active_dofs, const MPI_Comm mpi_communicator, const bool verbose=false) const
void resolve_indices(std::vector< types::global_dof_index > &indices) const
bool has_inhomogeneities() const
static const Table< 2, bool > default_empty_table
std::size_t memory_consumption() const
void add_lines(const std::vector< bool > &lines)
size_type max_constraint_indirections() const
size_type calculate_line_index(const size_type line) const
static ::ExceptionBase & ExcMatrixNotClosed()
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
static ::ExceptionBase & ExcNotQuadratic()
void swap(Vector< Number > &u, Vector< Number > &v)
void distribute(VectorType &vec) const
void write_dot(std::ostream &) const
void copy_from(const ConstraintMatrix &other)
std::atomic< unsigned int > counter
size_type row_length(const size_type row) const
bool is_identity_constrained(const size_type index) const
void shift(const size_type offset)
bool are_identity_constrained(const size_type index1, const size_type index2) const
static ::ExceptionBase & ExcNotImplemented()
bool is_element(const size_type index) const
void reinit(const IndexSet &local_constraints=IndexSet())
void add_selected_constraints(const ConstraintMatrix &constraints_in, const IndexSet &filter)
const std::vector< std::pair< size_type, double > > * get_constraint_entries(const size_type line) const
double get_inhomogeneity(const size_type line) const
void set_inhomogeneity(const size_type line, const double value)
const LineRange get_lines() const
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) const
bool is_compressed() const
std::map< unsigned int, T > some_to_some(const MPI_Comm &comm, const std::map< unsigned int, T > &objects_to_send)
const BlockIndices & get_column_indices() const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcMatrixIsClosed()
static ::ExceptionBase & ExcInternalError()