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_vector.h> 30 #include <deal.II/lac/petsc_block_vector.h> 31 #include <deal.II/lac/petsc_sparse_matrix.h> 32 #include <deal.II/lac/petsc_block_sparse_matrix.h> 33 #include <deal.II/lac/petsc_parallel_vector.h> 34 #include <deal.II/lac/petsc_parallel_block_vector.h> 35 #include <deal.II/lac/petsc_parallel_sparse_matrix.h> 36 #include <deal.II/lac/petsc_parallel_block_sparse_matrix.h> 37 #include <deal.II/lac/trilinos_vector.h> 38 #include <deal.II/lac/trilinos_block_vector.h> 39 #include <deal.II/lac/trilinos_sparse_matrix.h> 40 #include <deal.II/lac/trilinos_block_sparse_matrix.h> 41 #include <deal.II/lac/matrix_block.h> 42 #include <deal.II/lac/diagonal_matrix.h> 49 DEAL_II_NAMESPACE_OPEN
61 return (p.second == 0);
77 return line == a.
line;
95 for (std::set<size_type>::const_iterator
106 if (
lines[i] ==
true)
124 const std::vector<std::pair<size_type,double> > &col_val_pairs)
137 for (std::vector<std::pair<size_type,double> >::const_iterator
138 col_val_pair = col_val_pairs.begin();
139 col_val_pair!=col_val_pairs.end(); ++col_val_pair)
141 Assert (line != col_val_pair->first,
142 ExcMessage (
"Can't constrain a degree of freedom to itself"));
144 for (ConstraintLine::Entries::const_iterator
146 p != line_ptr->
entries.end(); ++p)
147 if (p->first == col_val_pair->first)
150 Assert (p->second == col_val_pair->second,
152 p->second, col_val_pair->second));
156 line_ptr->
entries.push_back (*col_val_pair);
170 ExcMessage (
"Filter needs to be larger than constraint matrix size."));
171 for (std::vector<ConstraintLine>::const_iterator line=constraints.
lines.begin();
172 line!=constraints.
lines.end(); ++line)
178 for (
size_type i=0; i<line->entries.size(); ++i)
179 if (filter.
is_element(line->entries[i].first))
181 line->entries[i].second);
198 std::vector<size_type> new_lines (
lines_cache.size(),
201 for (std::vector<ConstraintLine>::const_iterator line=
lines.begin();
214 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
215 line!=
lines.end(); ++line)
219 line->entries.erase (std::remove_if (line->entries.begin(),
222 line->entries.end());
235 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
236 line!=
lines.end(); ++line)
238 for (ConstraintLine::Entries::iterator it = line->entries.begin(); it!=line->entries.end(); ++it)
240 largest_idx=std::max(largest_idx, it->first);
260 bool chained_constraint_replaced =
false;
262 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
263 line!=
lines.end(); ++line)
276 while (entry < line->entries.size())
284 chained_constraint_replaced =
true;
287 const size_type dof_index = line->entries[entry].first;
288 const double weight = line->entries[entry].second;
290 Assert (dof_index != line->line,
291 ExcMessage (
"Cycle in constraints detected!"));
306 if (constrained_line->
entries.size() > 0)
310 ExcMessage (
"Cycle in constraints detected!"));
314 line->entries[entry] =
315 std::make_pair (constrained_line->
entries[0].first,
316 constrained_line->
entries[0].second *
321 .push_back (std::make_pair (constrained_line->
entries[i].first,
322 constrained_line->
entries[i].second *
330 Assert(n_replacements/2<largest_idx,
ExcMessage(
"Cycle in constraints detected!"));
331 if (n_replacements/2>=largest_idx)
342 line->entries.erase (line->entries.begin()+entry);
360 if (chained_constraint_replaced ==
false)
374 for (std::vector<ConstraintLine>::iterator line =
lines.begin();
375 line!=
lines.end(); ++line)
377 std::sort (line->entries.begin(), line->entries.end());
384 for (
size_type i=1; i<line->entries.size(); ++i)
385 if (line->entries[i].first == line->entries[i-1].first)
388 if (duplicates > 0 || line->entries.size() < line->entries.capacity())
395 new_entries = line->entries;
400 new_entries.reserve (line->entries.size() - duplicates);
401 new_entries.push_back(line->entries[0]);
402 for (
size_type j=1; j<line->entries.size(); ++j)
403 if (line->entries[j].first == line->entries[j-1].first)
405 Assert (new_entries.back().first == line->entries[j].first,
407 new_entries.back().second += line->entries[j].second;
410 new_entries.push_back (line->entries[j]);
412 Assert (new_entries.size() == line->entries.size() - duplicates,
417 for (
size_type j=1; j<new_entries.size(); ++j)
419 Assert (new_entries[j].first != new_entries[j-1].first,
421 Assert (new_entries[j].first > new_entries[j-1].first,
427 line->entries.swap (new_entries);
441 for (
size_type i=0; i<line->entries.size(); ++i)
442 sum += line->entries[i].second;
443 if ((sum != 1.0) && (std::fabs (sum-1.) < 1.e-13))
445 for (
size_type i=0; i<line->entries.size(); ++i)
446 line->entries[i].second /= sum;
447 line->inhomogeneity /= sum;
455 for (std::vector<ConstraintLine>::const_iterator line=
lines.begin();
456 line!=
lines.end(); ++line)
457 for (ConstraintLine::Entries::const_iterator
458 entry=line->entries.begin();
459 entry!=line->entries.end(); ++entry)
466 Assert (is_circle ==
false,
479 const bool allow_different_local_lines)
481 (void) allow_different_local_lines;
482 Assert(allow_different_local_lines ||
484 ExcMessage(
"local_lines for this and the other objects are not the same " 485 "although allow_different_local_lines is false."));
488 const bool object_was_sorted =
sorted;
500 for (std::vector<ConstraintLine>::iterator line=
lines.begin();
501 line!=
lines.end(); ++line)
504 for (
size_type i=0; i<line->entries.size(); ++i)
516 tmp.push_back(line->entries[i]);
527 const double weight = line->entries[i].second;
529 for (ConstraintLine::Entries::const_iterator j=other_line->begin();
530 j!=other_line->end(); ++j)
531 tmp.push_back (std::pair<size_type,double>(j->first,
540 line->entries.swap (tmp);
553 for (std::vector<ConstraintLine>::const_iterator line =
lines.begin();
554 line !=
lines.end(); ++line)
563 for (std::vector<ConstraintLine>::const_iterator line = other_constraints.
lines.begin();
564 line != other_constraints.
lines.end(); ++line)
570 lines.push_back(*line);
576 lines.push_back(*line);
583 switch (merge_conflict_behavior)
614 if (object_was_sorted ==
true)
633 for (std::vector<ConstraintLine>::iterator i =
lines.begin();
634 i !=
lines.end(); ++i)
637 for (ConstraintLine::Entries::iterator
638 j = i->entries.begin();
639 j != i->entries.end(); ++j)
658 std::vector<ConstraintLine> tmp;
663 std::vector<size_type> tmp;
719 ((entry != sparsity.
end(row)) &&
720 entry->is_valid_entry());
723 const size_type column = entry->column();
744 const size_type column = entry->column();
822 sparsity.
add (row, new_col);
825 if ((new_col < column) && (old_rowlength != new_rowlength))
827 old_rowlength = new_rowlength;
850 .entries.size(); ++q)
890 const std::pair<size_type,size_type>
891 block_index = index_mapping.global_to_local(row);
892 const size_type block_row = block_index.first;
901 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
904 block_sparsity = sparsity.
block(block_row, block_col);
907 entry = block_sparsity.
begin(block_index.second);
908 (entry != block_sparsity.
end(block_index.second)) &&
913 = index_mapping.local_to_global(block_col, entry->column());
931 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
934 block_sparsity = sparsity.
block(block_row,block_col);
937 entry = block_sparsity.
begin(block_index.second);
938 (entry != block_sparsity.
end(block_index.second)) &&
943 = index_mapping.local_to_global (block_col, entry->column());
1001 const std::pair<size_type,size_type>
1002 block_index = index_mapping.global_to_local(row);
1003 const size_type block_row = block_index.first;
1004 const size_type local_row = block_index.second;
1019 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1022 block_sparsity = sparsity.
block(block_row, block_col);
1027 = index_mapping.local_to_global(block_col,
1036 .entries.size(); ++q)
1047 for (
size_type block_col=0; block_col<n_blocks; ++block_col)
1050 block_sparsity = sparsity.
block(block_row,block_col);
1055 = index_mapping.local_to_global (block_col,
1095 return ((p.
entries.size() == 1) &&
1096 (p.
entries[0].second == 1.0));
1110 return ((p.
entries.size() == 1) &&
1111 (p.
entries[0].first == index2) &&
1112 (p.
entries[0].second == 1.0));
1121 return ((p.
entries.size() == 1) &&
1122 (p.
entries[0].first == index1) &&
1123 (p.
entries[0].second == 1.0));
1135 for (std::vector<ConstraintLine>::const_iterator i=
lines.begin();
1136 i!=
lines.end(); ++i)
1139 return_value = std::max(return_value,
1140 static_cast<size_type>(i->entries.size()));
1142 return return_value;
1149 for (std::vector<ConstraintLine>::const_iterator i=
lines.begin();
1150 i!=
lines.end(); ++i)
1151 if (i->inhomogeneity != 0.)
1163 if (
lines[i].entries.size() > 0)
1166 out <<
" " <<
lines[i].line
1167 <<
" " <<
lines[i].entries[j].first
1168 <<
": " <<
lines[i].entries[j].second <<
"\n";
1171 if (
lines[i].inhomogeneity != 0)
1172 out <<
" " <<
lines[i].line
1173 <<
": " <<
lines[i].inhomogeneity <<
"\n";
1180 if (
lines[i].inhomogeneity != 0)
1181 out <<
" " <<
lines[i].line
1182 <<
" = " <<
lines[i].inhomogeneity
1185 out <<
" " <<
lines[i].line <<
" = 0\n";
1197 out <<
"digraph constraints {" 1202 if (
lines[i].entries.size() > 0)
1204 out <<
" " <<
lines[i].line <<
"->" <<
lines[i].entries[j].first
1206 <<
lines[i].entries[j].second
1209 out <<
" " <<
lines[i].line <<
"\n";
1211 out <<
"}" << std::endl;
1230 const unsigned int indices_size = indices.size();
1231 const std::vector<std::pair<types::global_dof_index,double> > *line_ptr;
1232 for (
unsigned int i=0; i<indices_size; ++i)
1239 const unsigned int line_size = line_ptr->size();
1240 for (
unsigned int j=0; j<line_size; ++j)
1241 indices.push_back((*line_ptr)[j].first);
1246 std::sort(indices.begin(),indices.end());
1247 std::vector<types::global_dof_index>::iterator it;
1248 it = std::unique(indices.begin(),indices.end());
1249 indices.resize(it-indices.begin());
1264 #define VECTOR_FUNCTIONS(VectorType) \ 1265 template void ConstraintMatrix::condense<VectorType >(const VectorType &uncondensed,\ 1266 VectorType &condensed) const;\ 1267 template void ConstraintMatrix::condense<VectorType >(VectorType &vec) const;\ 1268 template void ConstraintMatrix:: \ 1269 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1270 const std::vector<ConstraintMatrix::size_type> &, \ 1272 const FullMatrix<VectorType::value_type> &) const;\ 1273 template void ConstraintMatrix:: \ 1274 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1275 const std::vector<ConstraintMatrix::size_type> &, \ 1276 const std::vector<ConstraintMatrix::size_type> &, \ 1278 const FullMatrix<VectorType::value_type> &, \ 1281 #define PARALLEL_VECTOR_FUNCTIONS(VectorType) \ 1282 template void ConstraintMatrix:: \ 1283 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1284 const std::vector<ConstraintMatrix::size_type> &, \ 1286 const FullMatrix<VectorType::value_type> &) const;\ 1287 template void ConstraintMatrix:: \ 1288 distribute_local_to_global<VectorType > (const Vector<VectorType::value_type> &, \ 1289 const std::vector<ConstraintMatrix::size_type> &, \ 1290 const std::vector<ConstraintMatrix::size_type> &, \ 1292 const FullMatrix<VectorType::value_type> &, \ 1295 #ifdef DEAL_II_WITH_PETSC 1300 #ifdef DEAL_II_WITH_TRILINOS 1305 #define MATRIX_VECTOR_FUNCTIONS(MatrixType, VectorType) \ 1306 template void ConstraintMatrix:: \ 1307 distribute_local_to_global<MatrixType,VectorType > (const FullMatrix<MatrixType::value_type> &, \ 1308 const Vector<VectorType::value_type> &, \ 1309 const std::vector<ConstraintMatrix::size_type> &, \ 1313 internal::bool2type<false>) const 1314 #define MATRIX_FUNCTIONS(MatrixType) \ 1315 template void ConstraintMatrix:: \ 1316 distribute_local_to_global<MatrixType,Vector<MatrixType::value_type> > (const FullMatrix<MatrixType::value_type> &, \ 1317 const Vector<MatrixType::value_type> &, \ 1318 const std::vector<ConstraintMatrix::size_type> &, \ 1320 Vector<MatrixType::value_type> &, \ 1322 internal::bool2type<false>) const 1323 #define BLOCK_MATRIX_VECTOR_FUNCTIONS(MatrixType, VectorType) \ 1324 template void ConstraintMatrix:: \ 1325 distribute_local_to_global<MatrixType,VectorType > (const FullMatrix<MatrixType::value_type> &, \ 1326 const Vector<VectorType::value_type> &, \ 1327 const std::vector<ConstraintMatrix::size_type> &, \ 1331 internal::bool2type<true>) const 1332 #define BLOCK_MATRIX_FUNCTIONS(MatrixType) \ 1333 template void ConstraintMatrix:: \ 1334 distribute_local_to_global<MatrixType,Vector<MatrixType::value_type> > (const FullMatrix<MatrixType::value_type> &, \ 1335 const Vector<MatrixType::value_type> &, \ 1336 const std::vector<ConstraintMatrix::size_type> &, \ 1338 Vector<MatrixType::value_type> &, \ 1340 internal::bool2type<true>) const 1346 MATRIX_FUNCTIONS(
FullMatrix<std::complex<double> >);
1363 #ifdef DEAL_II_WITH_PETSC 1374 #ifdef DEAL_II_WITH_TRILINOS 1384 #define SPARSITY_FUNCTIONS(SparsityPatternType) \ 1385 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1386 const std::vector<ConstraintMatrix::size_type> &, \ 1387 SparsityPatternType &, \ 1389 const Table<2,bool> &, \ 1390 internal::bool2type<false>) const; \ 1391 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1392 const std::vector<ConstraintMatrix::size_type> &, \ 1393 const std::vector<ConstraintMatrix::size_type> &, \ 1394 SparsityPatternType &, \ 1396 const Table<2,bool> &) const 1397 #define BLOCK_SPARSITY_FUNCTIONS(SparsityPatternType) \ 1398 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1399 const std::vector<ConstraintMatrix::size_type> &, \ 1400 SparsityPatternType &, \ 1402 const Table<2,bool> &, \ 1403 internal::bool2type<true>) const; \ 1404 template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \ 1405 const std::vector<ConstraintMatrix::size_type> &, \ 1406 const std::vector<ConstraintMatrix::size_type> &, \ 1407 SparsityPatternType &, \ 1409 const Table<2,bool> &) const 1416 #ifdef DEAL_II_WITH_TRILINOS 1422 #define ONLY_MATRIX_FUNCTIONS(MatrixType) \ 1423 template void ConstraintMatrix::distribute_local_to_global<MatrixType > (\ 1424 const FullMatrix<MatrixType::value_type> &, \ 1425 const std::vector<ConstraintMatrix::size_type> &, \ 1426 const std::vector<ConstraintMatrix::size_type> &, \ 1438 #ifdef DEAL_II_WITH_TRILINOS 1443 #ifdef DEAL_II_WITH_PETSC 1450 #include "constraint_matrix.inst" 1458 #define SCRATCH_INITIALIZER(Number,Name) \ 1459 ConstraintMatrixData<Number>::ScratchData scratch_data_initializer_##Name; \ 1460 template<> Threads::ThreadLocalStorage<ConstraintMatrixData<Number>::ScratchData> \ 1461 ConstraintMatrixData<Number>::scratch_data(scratch_data_initializer_##Name) 1463 SCRATCH_INITIALIZER(
double,
double);
1464 SCRATCH_INITIALIZER(
float,
float);
1465 SCRATCH_INITIALIZER(std::complex<double>,cdouble);
1466 SCRATCH_INITIALIZER(std::complex<float>,cfloat);
1467 #undef SCRATCH_INITIALIZER 1471 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)
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
DEAL_VOLATILE unsigned int counter
void add(const size_type i, const size_type j)
size_type n_block_cols() const
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)
#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
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 distribute(VectorType &vec) const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void write_dot(std::ostream &) const
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)
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) const
bool is_compressed() const
const BlockIndices & get_column_indices() const
static ::ExceptionBase & ExcMatrixIsClosed()
static ::ExceptionBase & ExcInternalError()