Reference documentation for deal.II version 9.0.0
constraint_matrix.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/constraint_matrix.h>
17 #include <deal.II/lac/constraint_matrix.templates.h>
18 
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>
40 
41 #include <algorithm>
42 #include <numeric>
43 #include <set>
44 #include <ostream>
45 #include <boost/serialization/utility.hpp>
46 
47 DEAL_II_NAMESPACE_OPEN
48 
49 
50 
51 // Static member variable
53 
54 
55 
56 void
58 {
59  lines = other.lines;
60  lines_cache = other.lines_cache;
61  local_lines = other.local_lines;
62  sorted = other.sorted;
63 }
64 
65 
66 
67 bool
68 ConstraintMatrix::check_zero_weight (const std::pair<size_type, double> &p)
69 {
70  return (p.second == 0);
71 }
72 
73 
74 
75 bool
77 {
78  return index < a.index;
79 }
80 
81 
82 
83 bool
85 {
86  return index == a.index;
87 }
88 
89 
90 
91 std::size_t
93 {
97 }
98 
99 
100 
102 {
103  return boost::make_iterator_range(lines.begin(), lines.end());
104 }
105 
106 
107 
108 bool ConstraintMatrix::is_consistent_in_parallel(const std::vector<IndexSet> &locally_owned_dofs,
109  const IndexSet &locally_active_dofs,
110  const MPI_Comm mpi_communicator,
111  const bool verbose) const
112 {
113  ConstraintLine empty;
114  empty.inhomogeneity = 0.0;
115 
116  // Helper to return a reference to the ConstraintLine object that belongs to row @p row.
117  // We don't want to make copies but to return a reference, we need an empty object that
118  // we store above.
119  auto get_line = [&] (const size_type row) -> const ConstraintLine&
120  {
121  const size_type line_index = calculate_line_index(row);
122  if (line_index >= lines_cache.size() ||
124  {
125  empty.index = row;
126  return empty;
127  }
128  else
129  return lines[lines_cache[line_index]];
130  };
131 
132  // identify non-owned rows and send to owner:
133  std::map< unsigned int, std::vector<ConstraintLine> > to_send;
134 
135  const unsigned int myid = ::Utilities::MPI::this_mpi_process(mpi_communicator);
136  const unsigned int nproc = ::Utilities::MPI::n_mpi_processes(mpi_communicator);
137 
138  // We will send all locally active dofs that are not locally owned for checking. Note
139  // that we allow constraints to differ on locally_relevant (and not active) DoFs.
140  IndexSet non_owned = locally_active_dofs;
141  non_owned.subtract_set(locally_owned_dofs[myid]);
142  for (unsigned int owner=0; owner<nproc; ++owner)
143  {
144  // find all lines to send to @p owner
145  IndexSet indices_to_send = non_owned & locally_owned_dofs[owner];
146  for (const auto &row_idx : indices_to_send)
147  {
148  to_send[owner].push_back(get_line(row_idx));
149  }
150  }
151 
152  std::map<unsigned int, std::vector<ConstraintLine> > received = Utilities::MPI::some_to_some (mpi_communicator, to_send);
153 
154  unsigned int inconsistent = 0;
155 
156  // from each processor:
157  for (const auto &kv : received)
158  {
159  // for each incoming line:
160  for (auto &lineit : kv.second)
161  {
162  const ConstraintLine &reference = get_line(lineit.index);
163 
164  if (lineit.inhomogeneity != reference.inhomogeneity)
165  {
166  ++inconsistent;
167 
168  if (verbose)
169  std::cout << "Proc " << myid
170  << " got line " << lineit.index
171  << " from " << kv.first
172  << " inhomogeneity " << lineit.inhomogeneity << " != " << reference.inhomogeneity << std::endl;
173  }
174  else if (lineit.entries != reference.entries)
175  {
176  ++inconsistent;
177  if (verbose)
178  std::cout << "Proc " << myid
179  << " got line " << lineit.index
180  << " from " << kv.first
181  << " wrong values!"
182  << std::endl;
183  }
184  }
185  }
186 
187  const unsigned int total = Utilities::MPI::sum(inconsistent, mpi_communicator);
188  if (verbose && total>0 && myid==0)
189  std::cout << total << " inconsistent lines discovered!" << std::endl;
190  return total==0;
191 }
192 
193 
194 
195 void
196 ConstraintMatrix::add_lines (const std::set<size_type> &lines)
197 {
198  for (std::set<size_type>::const_iterator
199  i = lines.begin(); i != lines.end(); ++i)
200  add_line (*i);
201 }
202 
203 
204 
205 void
206 ConstraintMatrix::add_lines (const std::vector<bool> &lines)
207 {
208  for (size_type i=0; i<lines.size(); ++i)
209  if (lines[i] == true)
210  add_line (i);
211 }
212 
213 
214 
215 void
217 {
218  for (size_type i=0; i<lines.n_elements(); ++i)
219  add_line (lines.nth_index_in_set(i));
220 }
221 
222 
223 
224 void
226 (const size_type line,
227  const std::vector<std::pair<size_type,double> > &col_val_pairs)
228 {
229  Assert (sorted==false, ExcMatrixIsClosed());
230  Assert (is_constrained(line), ExcLineInexistant(line));
231 
233  Assert (line_ptr->index == line, ExcInternalError());
234 
235  // if in debug mode, check whether an entry for this column already
236  // exists and if its the same as the one entered at present
237  //
238  // in any case: skip this entry if an entry for this column already
239  // exists, since we don't want to enter it twice
240  for (std::vector<std::pair<size_type,double> >::const_iterator
241  col_val_pair = col_val_pairs.begin();
242  col_val_pair!=col_val_pairs.end(); ++col_val_pair)
243  {
244  Assert (line != col_val_pair->first,
245  ExcMessage ("Can't constrain a degree of freedom to itself"));
246 
247  for (ConstraintLine::Entries::const_iterator
248  p=line_ptr->entries.begin();
249  p != line_ptr->entries.end(); ++p)
250  if (p->first == col_val_pair->first)
251  {
252  // entry exists, break innermost loop
253  Assert (p->second == col_val_pair->second,
254  ExcEntryAlreadyExists(line, col_val_pair->first,
255  p->second, col_val_pair->second));
256  break;
257  }
258 
259  line_ptr->entries.push_back (*col_val_pair);
260  }
261 }
262 
263 
264 
266 (const ConstraintMatrix &constraints,
267  const IndexSet &filter)
268 {
269  if (constraints.n_constraints() == 0)
270  return;
271 
272  Assert (filter.size() > constraints.lines.back().index,
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)
276  if (filter.is_element(line->index))
277  {
278  const size_type row = filter.index_within_set (line->index);
279  add_line (row);
280  set_inhomogeneity (row, line->inhomogeneity);
281  for (size_type i=0; i<line->entries.size(); ++i)
282  if (filter.is_element(line->entries[i].first))
283  add_entry (row, filter.index_within_set (line->entries[i].first),
284  line->entries[i].second);
285  }
286 }
287 
288 
289 
291 {
292  if (sorted == true)
293  return;
294 
295  // sort the lines
296  std::sort (lines.begin(), lines.end());
297 
298  // update list of pointers and give the vector a sharp size since we
299  // won't modify the size any more after this point.
300  {
301  std::vector<size_type> new_lines (lines_cache.size(),
303  size_type counter = 0;
304  for (std::vector<ConstraintLine>::const_iterator line=lines.begin();
305  line!=lines.end(); ++line, ++counter)
306  new_lines[calculate_line_index(line->index)] = counter;
307  std::swap (lines_cache, new_lines);
308  }
309 
310  // in debug mode: check whether we really set the pointers correctly.
311  for (size_type i=0; i<lines_cache.size(); ++i)
313  Assert (i == calculate_line_index(lines[lines_cache[i]].index),
314  ExcInternalError());
315 
316  // first, strip zero entries, as we have to do that only once
317  for (std::vector<ConstraintLine>::iterator line = lines.begin();
318  line!=lines.end(); ++line)
319  // first remove zero entries. that would mean that in the linear
320  // constraint for a node, x_i = ax_1 + bx_2 + ..., another node times 0
321  // appears. obviously, 0*something can be omitted
322  line->entries.erase (std::remove_if (line->entries.begin(),
323  line->entries.end(),
325  line->entries.end());
326 
327 
328 
329 #ifdef DEBUG
330  // In debug mode we are computing an estimate for the maximum number
331  // of constraints so that we can bail out if there is a cycle in the
332  // constraints (which is easier than searching for cycles in the graph).
333  //
334  // Let us figure out the largest dof index. This is an upper bound for the
335  // number of constraints because it is an approximation for the number of dofs
336  // in our system.
337  size_type largest_idx = 0;
338  for (std::vector<ConstraintLine>::iterator line = lines.begin();
339  line!=lines.end(); ++line)
340  {
341  for (ConstraintLine::Entries::iterator it = line->entries.begin(); it!=line->entries.end(); ++it)
342  {
343  largest_idx=std::max(largest_idx, it->first);
344  }
345  }
346 #endif
347 
348  // replace references to dofs that are themselves constrained. note that
349  // because we may replace references to other dofs that may themselves be
350  // constrained to third ones, we have to iterate over all this until we
351  // replace no chains of constraints any more
352  //
353  // the iteration replaces references to constrained degrees of freedom by
354  // second-order references. for example if x3=x0/2+x2/2 and x2=x0/2+x1/2,
355  // then the new list will be x3=x0/2+x0/4+x1/4. note that x0 appear
356  // twice. we will throw this duplicate out in the following step, where
357  // we sort the list so that throwing out duplicates becomes much more
358  // efficient. also, we have to do it only once, rather than in each
359  // iteration
360  size_type iteration = 0;
361  while (true)
362  {
363  bool chained_constraint_replaced = false;
364 
365  for (std::vector<ConstraintLine>::iterator line = lines.begin();
366  line!=lines.end(); ++line)
367  {
368 #ifdef DEBUG
369  // we need to keep track of how many replacements we do in this line, because we can
370  // end up in a cycle A->B->C->A without the number of entries growing.
371  size_type n_replacements = 0;
372 #endif
373 
374  // loop over all entries of this line (including ones that we
375  // have appended in this go around) and see whether they are
376  // further constrained. ignore elements that we don't store on
377  // the current processor
378  size_type entry = 0;
379  while (entry < line->entries.size())
380  if (((local_lines.size() == 0)
381  ||
382  (local_lines.is_element(line->entries[entry].first)))
383  &&
384  is_constrained (line->entries[entry].first))
385  {
386  // ok, this entry is further constrained:
387  chained_constraint_replaced = true;
388 
389  // look up the chain of constraints for this entry
390  const size_type dof_index = line->entries[entry].first;
391  const double weight = line->entries[entry].second;
392 
393  Assert (dof_index != line->index,
394  ExcMessage ("Cycle in constraints detected!"));
395 
396  const ConstraintLine *constrained_line =
397  &lines[lines_cache[calculate_line_index(dof_index)]];
398  Assert (constrained_line->index == dof_index,
399  ExcInternalError());
400 
401  // now we have to replace an entry by its expansion. we do
402  // that by overwriting the entry by the first entry of the
403  // expansion and adding the remaining ones to the end,
404  // where we will later process them once more
405  //
406  // we can of course only do that if the DoF that we are
407  // currently handle is constrained by a linear combination
408  // of other dofs:
409  if (constrained_line->entries.size() > 0)
410  {
411  for (size_type i=0; i<constrained_line->entries.size(); ++i)
412  Assert (dof_index != constrained_line->entries[i].first,
413  ExcMessage ("Cycle in constraints detected!"));
414 
415  // replace first entry, then tack the rest to the end
416  // of the list
417  line->entries[entry] =
418  std::make_pair (constrained_line->entries[0].first,
419  constrained_line->entries[0].second *
420  weight);
421 
422  for (size_type i=1; i<constrained_line->entries.size(); ++i)
423  line->entries.emplace_back (constrained_line->entries[i].first,
424  constrained_line->entries[i].second
425  * weight);
426 
427 #ifdef DEBUG
428  // keep track of how many entries we replace in this
429  // line. If we do more than there are constraints or
430  // dofs in our system, we must have a cycle.
431  ++n_replacements;
432  Assert(n_replacements/2<largest_idx, ExcMessage("Cycle in constraints detected!"));
433  if (n_replacements/2>=largest_idx)
434  return; // this enables us to test for this Exception.
435 #endif
436  }
437  else
438  // the DoF that we encountered is not constrained by a
439  // linear combination of other dofs but is equal to just
440  // the inhomogeneity (i.e. its chain of entries is
441  // empty). in that case, we can't just overwrite the
442  // current entry, but we have to actually eliminate it
443  {
444  line->entries.erase (line->entries.begin()+entry);
445  }
446 
447  line->inhomogeneity += constrained_line->inhomogeneity *
448  weight;
449 
450  // now that we're here, do not increase index by one but
451  // rather make another pass for the present entry because
452  // we have replaced the present entry by another one, or
453  // because we have deleted it and shifted all following
454  // ones one forward
455  }
456  else
457  // entry not further constrained. just move ahead by one
458  ++entry;
459  }
460 
461  // if we didn't do anything in this round, then quit the loop
462  if (chained_constraint_replaced == false)
463  break;
464 
465  // increase iteration count. note that we should not iterate more
466  // times than there are constraints, since this puts a natural upper
467  // bound on the length of constraint chains
468  ++iteration;
469  Assert (iteration <= lines.size(), ExcInternalError());
470  }
471 
472  // finally sort the entries and re-scale them if necessary. in this step,
473  // we also throw out duplicates as mentioned above. moreover, as some
474  // entries might have had zero weights, we replace them by a vector with
475  // sharp sizes.
476  for (std::vector<ConstraintLine>::iterator line = lines.begin();
477  line!=lines.end(); ++line)
478  {
479  std::sort (line->entries.begin(), line->entries.end());
480 
481  // loop over the now sorted list and see whether any of the entries
482  // references the same dofs more than once in order to find how many
483  // non-duplicate entries we have. This lets us allocate the correct
484  // amount of memory for the constraint entries.
485  size_type duplicates = 0;
486  for (size_type i=1; i<line->entries.size(); ++i)
487  if (line->entries[i].first == line->entries[i-1].first)
488  duplicates++;
489 
490  if (duplicates > 0 || line->entries.size() < line->entries.capacity())
491  {
492  ConstraintLine::Entries new_entries;
493 
494  // if we have no duplicates, copy verbatim the entries. this way,
495  // the final size is of the vector is correct.
496  if (duplicates == 0)
497  new_entries = line->entries;
498  else
499  {
500  // otherwise, we need to go through the list by and and
501  // resolve the duplicates
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)
506  {
507  Assert (new_entries.back().first == line->entries[j].first,
508  ExcInternalError());
509  new_entries.back().second += line->entries[j].second;
510  }
511  else
512  new_entries.push_back (line->entries[j]);
513 
514  Assert (new_entries.size() == line->entries.size() - duplicates,
515  ExcInternalError());
516 
517  // make sure there are really no duplicates left and that the
518  // list is still sorted
519  for (size_type j=1; j<new_entries.size(); ++j)
520  {
521  Assert (new_entries[j].first != new_entries[j-1].first,
522  ExcInternalError());
523  Assert (new_entries[j].first > new_entries[j-1].first,
524  ExcInternalError());
525  }
526  }
527 
528  // replace old list of constraints for this dof by the new one
529  line->entries.swap (new_entries);
530  }
531 
532  // finally do the following check: if the sum of weights for the
533  // constraints is close to one, but not exactly one, then rescale all
534  // the weights so that they sum up to 1. this adds a little numerical
535  // stability and avoids all sorts of problems where the actual value
536  // is close to, but not quite what we expected
537  //
538  // the case where the weights don't quite sum up happens when we
539  // compute the interpolation weights "on the fly", i.e. not from
540  // precomputed tables. in this case, the interpolation weights are
541  // also subject to round-off
542  double sum = 0;
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))
546  {
547  for (size_type i=0; i<line->entries.size(); ++i)
548  line->entries[i].second /= sum;
549  line->inhomogeneity /= sum;
550  }
551  } // end of loop over all constraint lines
552 
553 #ifdef DEBUG
554  // if in debug mode: check that no dof is constrained to another dof that
555  // is also constrained. exclude dofs from this check whose constraint
556  // lines are not stored on the local processor
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)
562  if ((local_lines.size() == 0)
563  ||
564  (local_lines.is_element(entry->first)))
565  {
566  // make sure that entry->first is not the index of a line itself
567  const bool is_circle = is_constrained(entry->first);
568  Assert (is_circle == false,
569  ExcDoFConstrainedToConstrainedDoF(line->index, entry->first));
570  }
571 #endif
572 
573  sorted = true;
574 }
575 
576 
577 
578 void
579 ConstraintMatrix::merge (const ConstraintMatrix &other_constraints,
580  const MergeConflictBehavior merge_conflict_behavior,
581  const bool allow_different_local_lines)
582 {
583  (void) allow_different_local_lines;
584  Assert(allow_different_local_lines ||
585  local_lines == other_constraints.local_lines,
586  ExcMessage("local_lines for this and the other objects are not the same "
587  "although allow_different_local_lines is false."));
588 
589  // store the previous state with respect to sorting
590  const bool object_was_sorted = sorted;
591  sorted = false;
592 
593  // first action is to fold into the present object possible constraints
594  // in the second object. we don't strictly need to do this any more since
595  // the ConstraintMatrix has learned to deal with chains of constraints in
596  // the close() function, but we have traditionally done this and it's not
597  // overly hard to do.
598  //
599  // for this, loop over all constraints and replace the constraint lines
600  // with a new one where constraints are replaced if necessary.
602  for (std::vector<ConstraintLine>::iterator line=lines.begin();
603  line!=lines.end(); ++line)
604  {
605  tmp.clear ();
606  for (size_type i=0; i<line->entries.size(); ++i)
607  {
608  // if the present dof is not stored, or not constrained, or if we won't take the
609  // constraint from the other object, then simply copy it over
610  if ((other_constraints.local_lines.size() != 0
611  && other_constraints.local_lines.is_element(line->entries[i].first) == false)
612  ||
613  other_constraints.is_constrained(line->entries[i].first) == false
614  ||
615  ((merge_conflict_behavior != right_object_wins)
616  && other_constraints.is_constrained(line->entries[i].first)
617  && this->is_constrained(line->entries[i].first)))
618  tmp.push_back(line->entries[i]);
619  else
620  // otherwise resolve further constraints by replacing the old
621  // entry by a sequence of new entries taken from the other
622  // object, but with multiplied weights
623  {
624  const ConstraintLine::Entries *other_line
625  = other_constraints.get_constraint_entries (line->entries[i].first);
626  Assert (other_line != nullptr,
627  ExcInternalError());
628 
629  const double weight = line->entries[i].second;
630 
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);
634 
635  line->inhomogeneity
636  += other_constraints.get_inhomogeneity(line->entries[i].first) *
637  weight;
638  }
639  }
640  // finally exchange old and newly resolved line
641  line->entries.swap (tmp);
642  }
643 
644  if (local_lines.size() != 0)
645  local_lines.add_indices(other_constraints.local_lines);
646 
647  {
648  // do not bother to resize the lines cache exactly since it is pretty
649  // cheap to adjust it along the way.
650  std::fill(lines_cache.begin(), lines_cache.end(), numbers::invalid_size_type);
651 
652  // reset lines_cache for our own constraints
653  size_type index = 0;
654  for (std::vector<ConstraintLine>::const_iterator line = lines.begin();
655  line != lines.end(); ++line)
656  {
657  size_type local_line_no = calculate_line_index(line->index);
658  if (local_line_no >= lines_cache.size())
659  lines_cache.resize(local_line_no+1, numbers::invalid_size_type);
660  lines_cache[local_line_no] = index++;
661  }
662 
663  // Add other_constraints to lines cache and our list of constraints
664  for (std::vector<ConstraintLine>::const_iterator line = other_constraints.lines.begin();
665  line != other_constraints.lines.end(); ++line)
666  {
667  const size_type local_line_no = calculate_line_index(line->index);
668  if (local_line_no >= lines_cache.size())
669  {
670  lines_cache.resize(local_line_no+1, numbers::invalid_size_type);
671  lines.push_back(*line);
672  lines_cache[local_line_no] = index++;
673  }
674  else if (lines_cache[local_line_no] == numbers::invalid_size_type)
675  {
676  // there are no constraints for that line yet
677  lines.push_back(*line);
678  AssertIndexRange(local_line_no, lines_cache.size());
679  lines_cache[local_line_no] = index++;
680  }
681  else
682  {
683  // we already store that line
684  switch (merge_conflict_behavior)
685  {
687  AssertThrow (false,
688  ExcDoFIsConstrainedFromBothObjects (line->index));
689  break;
690 
691  case left_object_wins:
692  // ignore this constraint
693  break;
694 
695  case right_object_wins:
696  AssertIndexRange(local_line_no, lines_cache.size());
697  lines[lines_cache[local_line_no]] = *line;
698  break;
699 
700  default:
701  Assert (false, ExcNotImplemented());
702  }
703  }
704  }
705 
706  // check that we set the pointers correctly
707  for (size_type i=0; i<lines_cache.size(); ++i)
709  Assert (i == calculate_line_index(lines[lines_cache[i]].index),
710  ExcInternalError());
711  }
712 
713  // if the object was sorted before, then make sure it is so afterward as
714  // well. otherwise leave everything in the unsorted state
715  if (object_was_sorted == true)
716  close ();
717 }
718 
719 
720 
722 {
723  if (local_lines.size() == 0)
724  lines_cache.insert (lines_cache.begin(), offset,
726  else
727  {
728  // shift local_lines
729  IndexSet new_local_lines(local_lines.size());
730  new_local_lines.add_indices(local_lines, offset);
731  std::swap(local_lines, new_local_lines);
732  }
733 
734  for (std::vector<ConstraintLine>::iterator i = lines.begin();
735  i != lines.end(); ++i)
736  {
737  i->index += offset;
738  for (ConstraintLine::Entries::iterator
739  j = i->entries.begin();
740  j != i->entries.end(); ++j)
741  j->first += offset;
742  }
743 
744 #ifdef DEBUG
745  // make sure that lines, lines_cache and local_lines
746  // are still linked correctly
747  for (size_type i=0; i<lines_cache.size(); ++i)
749  calculate_line_index(lines[lines_cache[i]].index) == i,
750  ExcInternalError());
751 #endif
752 }
753 
754 
755 
757 {
758  {
759  std::vector<ConstraintLine> tmp;
760  lines.swap (tmp);
761  }
762 
763  {
764  std::vector<size_type> tmp;
765  lines_cache.swap (tmp);
766  }
767 
768  sorted = false;
769 }
770 
771 
772 
773 void ConstraintMatrix::reinit (const IndexSet &local_constraints)
774 {
775  local_lines = local_constraints;
776 
777  // make sure the IndexSet is compressed. Otherwise this can lead to crashes
778  // that are hard to find (only happen in release mode).
779  // see tests/mpi/constraint_matrix_crash_01
781 
782  clear();
783 }
784 
785 
786 
788 {
789  Assert (sorted == true, ExcMatrixNotClosed());
790  Assert (sparsity.is_compressed() == false, ExcMatrixIsClosed());
791  Assert (sparsity.n_rows() == sparsity.n_cols(), ExcNotQuadratic());
792 
793  // store for each index whether it must be distributed or not. If entry
794  // is numbers::invalid_unsigned_int, no distribution is necessary.
795  // otherwise, the number states which line in the constraint matrix
796  // handles this index
797  std::vector<size_type> distribute(sparsity.n_rows(),
799 
800  for (size_type c=0; c<lines.size(); ++c)
801  distribute[lines[c].index] = c;
802 
803  const size_type n_rows = sparsity.n_rows();
804  for (size_type row=0; row<n_rows; ++row)
805  {
807  {
808  // regular line. loop over cols all valid cols. note that this
809  // changes the line we are presently working on: we add additional
810  // entries. these are put to the end of the row. however, as
811  // constrained nodes cannot be constrained to other constrained
812  // nodes, nothing will happen if we run into these added nodes, as
813  // they can't be distributed further. we might store the position of
814  // the last old entry and stop work there, but since operating on
815  // the newly added ones only takes two comparisons (column index
816  // valid, distribute[column] necessarily
817  // ==numbers::invalid_size_type), it is cheaper to not do so and
818  // run right until the end of the line
819  for (SparsityPattern::iterator entry = sparsity.begin(row);
820  ((entry != sparsity.end(row)) &&
821  entry->is_valid_entry());
822  ++entry)
823  {
824  const size_type column = entry->column();
825 
826  if (distribute[column] != numbers::invalid_size_type)
827  {
828  // distribute entry at regular row @p{row} and irregular
829  // column sparsity.colnums[j]
830  for (size_type q=0;
831  q!=lines[distribute[column]].entries.size();
832  ++q)
833  sparsity.add (row,
834  lines[distribute[column]].entries[q].first);
835  }
836  }
837  }
838  else
839  // row must be distributed. note that here the present row is not
840  // touched (unlike above)
841  {
842  for (SparsityPattern::iterator entry = sparsity.begin(row);
843  (entry != sparsity.end(row)) && entry->is_valid_entry(); ++entry)
844  {
845  const size_type column = entry->column();
846  if (distribute[column] == numbers::invalid_size_type)
847  // distribute entry at irregular row @p{row} and regular
848  // column sparsity.colnums[j]
849  for (size_type q=0;
850  q!=lines[distribute[row]].entries.size(); ++q)
851  sparsity.add (lines[distribute[row]].entries[q].first,
852  column);
853  else
854  // distribute entry at irregular row @p{row} and irregular
855  // column sparsity.get_column_numbers()[j]
856  for (size_type p=0; p!=lines[distribute[row]].entries.size(); ++p)
857  for (size_type q=0;
858  q!=lines[distribute[column]].entries.size(); ++q)
859  sparsity.add (lines[distribute[row]].entries[p].first,
860  lines[distribute[column]].entries[q].first);
861  }
862  }
863  }
864 
865  sparsity.compress();
866 }
867 
868 
869 
870 
872 {
873  Assert (sorted == true, ExcMatrixNotClosed());
874  Assert (sparsity.n_rows() == sparsity.n_cols(),
875  ExcNotQuadratic());
876 
877  // store for each index whether it must be distributed or not. If entry
878  // is numbers::invalid_unsigned_int, no distribution is necessary.
879  // otherwise, the number states which line in the constraint matrix
880  // handles this index
881  std::vector<size_type> distribute(sparsity.n_rows(),
883 
884  for (size_type c=0; c<lines.size(); ++c)
885  distribute[lines[c].index] = c;
886 
887  const size_type n_rows = sparsity.n_rows();
888  for (size_type row=0; row<n_rows; ++row)
889  {
891  // regular line. loop over cols. note that as we proceed to
892  // distribute cols, the loop may get longer
893  for (size_type j=0; j<sparsity.row_length(row); ++j)
894  {
895  const size_type column = sparsity.column_number(row,j);
896 
897  if (distribute[column] != numbers::invalid_size_type)
898  {
899  // distribute entry at regular row @p{row} and irregular
900  // column column. note that this changes the line we are
901  // presently working on: we add additional entries. if we
902  // add another entry at a column behind the present one, we
903  // will encounter it later on (but since it can't be
904  // further constrained, won't have to do anything about
905  // it). if we add it up front of the present column, we
906  // will find the present column later on again as it was
907  // shifted back (again nothing happens, in particular no
908  // endless loop, as when we encounter it the second time we
909  // won't be able to add more entries as they all already
910  // exist, but we do the same work more often than
911  // necessary, and the loop gets longer), so move the cursor
912  // one to the right in the case that we add an entry up
913  // front that did not exist before. check whether it
914  // existed before by tracking the length of this row
915  size_type old_rowlength = sparsity.row_length(row);
916  for (size_type q=0;
917  q!=lines[distribute[column]].entries.size();
918  ++q)
919  {
920  const size_type
921  new_col = lines[distribute[column]].entries[q].first;
922 
923  sparsity.add (row, new_col);
924 
925  const size_type new_rowlength = sparsity.row_length(row);
926  if ((new_col < column) && (old_rowlength != new_rowlength))
927  ++j;
928  old_rowlength = new_rowlength;
929  }
930  }
931  }
932  else
933  // row must be distributed
934  for (size_type j=0; j<sparsity.row_length(row); ++j)
935  {
936  const size_type column = sparsity.column_number(row,j);
937 
938  if (distribute[column] == numbers::invalid_size_type)
939  // distribute entry at irregular row @p{row} and regular
940  // column sparsity.colnums[j]
941  for (size_type q=0;
942  q!=lines[distribute[row]].entries.size(); ++q)
943  sparsity.add (lines[distribute[row]].entries[q].first,
944  column);
945  else
946  // distribute entry at irregular row @p{row} and irregular
947  // column sparsity.get_column_numbers()[j]
948  for (size_type p=0; p!=lines[distribute[row]].entries.size(); ++p)
949  for (size_type q=0;
950  q!=lines[distribute[sparsity.column_number(row,j)]]
951  .entries.size(); ++q)
952  sparsity.add (lines[distribute[row]].entries[p].first,
953  lines[distribute[sparsity.column_number(row,j)]]
954  .entries[q].first);
955  }
956  }
957 }
958 
959 
960 
962 {
963  Assert (sorted == true, ExcMatrixNotClosed());
964  Assert (sparsity.is_compressed() == false, ExcMatrixIsClosed());
965  Assert (sparsity.n_rows() == sparsity.n_cols(),
966  ExcNotQuadratic());
967  Assert (sparsity.n_block_rows() == sparsity.n_block_cols(),
968  ExcNotQuadratic());
969  Assert (sparsity.get_column_indices() == sparsity.get_row_indices(),
970  ExcNotQuadratic());
971 
972  const BlockIndices &
973  index_mapping = sparsity.get_column_indices();
974 
975  const size_type n_blocks = sparsity.n_block_rows();
976 
977  // store for each index whether it must be distributed or not. If entry
978  // is numbers::invalid_unsigned_int, no distribution is necessary.
979  // otherwise, the number states which line in the constraint matrix
980  // handles this index
981  std::vector<size_type> distribute (sparsity.n_rows(),
983 
984  for (size_type c=0; c<lines.size(); ++c)
985  distribute[lines[c].index] = c;
986 
987  const size_type n_rows = sparsity.n_rows();
988  for (size_type row=0; row<n_rows; ++row)
989  {
990  // get index of this row within the blocks
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;
994 
996  // regular line. loop over all columns and see whether this column
997  // must be distributed
998  {
999 
1000  // to loop over all entries in this row, we have to loop over all
1001  // blocks in this blockrow and the corresponding row therein
1002  for (size_type block_col=0; block_col<n_blocks; ++block_col)
1003  {
1004  const SparsityPattern &
1005  block_sparsity = sparsity.block(block_row, block_col);
1006 
1008  entry = block_sparsity.begin(block_index.second);
1009  (entry != block_sparsity.end(block_index.second)) &&
1010  entry->is_valid_entry();
1011  ++entry)
1012  {
1013  const size_type global_col
1014  = index_mapping.local_to_global(block_col, entry->column());
1015 
1016  if (distribute[global_col] != numbers::invalid_size_type)
1017  // distribute entry at regular row @p{row} and
1018  // irregular column global_col
1019  {
1020  for (size_type q=0;
1021  q!=lines[distribute[global_col]].entries.size(); ++q)
1022  sparsity.add (row,
1023  lines[distribute[global_col]].entries[q].first);
1024  }
1025  }
1026  }
1027  }
1028  else
1029  {
1030  // row must be distributed. split the whole row into the chunks
1031  // defined by the blocks
1032  for (size_type block_col=0; block_col<n_blocks; ++block_col)
1033  {
1034  const SparsityPattern &
1035  block_sparsity = sparsity.block(block_row,block_col);
1036 
1038  entry = block_sparsity.begin(block_index.second);
1039  (entry != block_sparsity.end(block_index.second)) &&
1040  entry->is_valid_entry();
1041  ++entry)
1042  {
1043  const size_type global_col
1044  = index_mapping.local_to_global (block_col, entry->column());
1045 
1046  if (distribute[global_col] == numbers::invalid_size_type)
1047  // distribute entry at irregular row @p{row} and
1048  // regular column global_col.
1049  {
1050  for (size_type q=0; q!=lines[distribute[row]].entries.size(); ++q)
1051  sparsity.add (lines[distribute[row]].entries[q].first, global_col);
1052  }
1053  else
1054  // distribute entry at irregular row @p{row} and
1055  // irregular column @p{global_col}
1056  {
1057  for (size_type p=0; p!=lines[distribute[row]].entries.size(); ++p)
1058  for (size_type q=0; q!=lines[distribute[global_col]].entries.size(); ++q)
1059  sparsity.add (lines[distribute[row]].entries[p].first,
1060  lines[distribute[global_col]].entries[q].first);
1061  }
1062  }
1063  }
1064  }
1065  }
1066 
1067  sparsity.compress();
1068 }
1069 
1070 
1071 
1072 
1074 {
1075  Assert (sorted == true, ExcMatrixNotClosed());
1076  Assert (sparsity.n_rows() == sparsity.n_cols(),
1077  ExcNotQuadratic());
1078  Assert (sparsity.n_block_rows() == sparsity.n_block_cols(),
1079  ExcNotQuadratic());
1080  Assert (sparsity.get_column_indices() == sparsity.get_row_indices(),
1081  ExcNotQuadratic());
1082 
1083  const BlockIndices &
1084  index_mapping = sparsity.get_column_indices();
1085 
1086  const size_type n_blocks = sparsity.n_block_rows();
1087 
1088  // store for each index whether it must be distributed or not. If entry
1089  // is numbers::invalid_unsigned_int, no distribution is necessary.
1090  // otherwise, the number states which line in the constraint matrix
1091  // handles this index
1092  std::vector<size_type> distribute (sparsity.n_rows(),
1094 
1095  for (size_type c=0; c<lines.size(); ++c)
1096  distribute[lines[c].index] = static_cast<signed int>(c);
1097 
1098  const size_type n_rows = sparsity.n_rows();
1099  for (size_type row=0; row<n_rows; ++row)
1100  {
1101  // get index of this row within the blocks
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;
1106 
1108  // regular line. loop over all columns and see whether this column
1109  // must be distributed. note that as we proceed to distribute cols,
1110  // the loop over cols may get longer.
1111  //
1112  // don't try to be clever here as in the algorithm for the
1113  // DynamicSparsityPattern, as that would be much more
1114  // complicated here. after all, we know that compressed patterns
1115  // are inefficient...
1116  {
1117 
1118  // to loop over all entries in this row, we have to loop over all
1119  // blocks in this blockrow and the corresponding row therein
1120  for (size_type block_col=0; block_col<n_blocks; ++block_col)
1121  {
1122  const DynamicSparsityPattern &
1123  block_sparsity = sparsity.block(block_row, block_col);
1124 
1125  for (size_type j=0; j<block_sparsity.row_length(local_row); ++j)
1126  {
1127  const size_type global_col
1128  = index_mapping.local_to_global(block_col,
1129  block_sparsity.column_number(local_row,j));
1130 
1131  if (distribute[global_col] != numbers::invalid_size_type)
1132  // distribute entry at regular row @p{row} and
1133  // irregular column global_col
1134  {
1135  for (size_type q=0;
1136  q!=lines[distribute[global_col]]
1137  .entries.size(); ++q)
1138  sparsity.add (row,
1139  lines[distribute[global_col]].entries[q].first);
1140  }
1141  }
1142  }
1143  }
1144  else
1145  {
1146  // row must be distributed. split the whole row into the chunks
1147  // defined by the blocks
1148  for (size_type block_col=0; block_col<n_blocks; ++block_col)
1149  {
1150  const DynamicSparsityPattern &
1151  block_sparsity = sparsity.block(block_row,block_col);
1152 
1153  for (size_type j=0; j<block_sparsity.row_length(local_row); ++j)
1154  {
1155  const size_type global_col
1156  = index_mapping.local_to_global (block_col,
1157  block_sparsity.column_number(local_row,j));
1158 
1159  if (distribute[global_col] == numbers::invalid_size_type)
1160  // distribute entry at irregular row @p{row} and
1161  // regular column global_col.
1162  {
1163  for (size_type q=0;
1164  q!=lines[distribute[row]].entries.size(); ++q)
1165  sparsity.add (lines[distribute[row]].entries[q].first,
1166  global_col);
1167  }
1168  else
1169  // distribute entry at irregular row @p{row} and
1170  // irregular column @p{global_col}
1171  {
1172  for (size_type p=0;
1173  p!=lines[distribute[row]].entries.size(); ++p)
1174  for (size_type q=0; q!=lines[distribute[global_col]].entries.size(); ++q)
1175  sparsity.add (lines[distribute[row]].entries[p].first,
1176  lines[distribute[global_col]].entries[q].first);
1177  }
1178  }
1179  }
1180  }
1181  }
1182 }
1183 
1184 
1185 
1187 {
1188  if (is_constrained(index) == false)
1189  return false;
1190 
1192  Assert (p.index == index, ExcInternalError());
1193 
1194  // return if an entry for this line was found and if it has only one
1195  // entry equal to 1.0
1196  return ((p.entries.size() == 1) &&
1197  (p.entries[0].second == 1.0));
1198 }
1199 
1200 
1202  const size_type index2) const
1203 {
1204  if (is_constrained(index1) == true)
1205  {
1207  Assert (p.index == index1, ExcInternalError());
1208 
1209  // return if an entry for this line was found and if it has only one
1210  // entry equal to 1.0 and that one is index2
1211  return ((p.entries.size() == 1) &&
1212  (p.entries[0].first == index2) &&
1213  (p.entries[0].second == 1.0));
1214  }
1215  else if (is_constrained(index2) == true)
1216  {
1218  Assert (p.index == index2, ExcInternalError());
1219 
1220  // return if an entry for this line was found and if it has only one
1221  // entry equal to 1.0 and that one is index1
1222  return ((p.entries.size() == 1) &&
1223  (p.entries[0].first == index1) &&
1224  (p.entries[0].second == 1.0));
1225  }
1226  else
1227  return false;
1228 }
1229 
1230 
1231 
1234 {
1235  size_type return_value = 0;
1236  for (std::vector<ConstraintLine>::const_iterator i=lines.begin();
1237  i!=lines.end(); ++i)
1238  // use static cast, since typeof(size)==std::size_t, which is !=
1239  // size_type on AIX
1240  return_value = std::max(return_value,
1241  static_cast<size_type>(i->entries.size()));
1242 
1243  return return_value;
1244 }
1245 
1246 
1247 
1249 {
1250  for (std::vector<ConstraintLine>::const_iterator i=lines.begin();
1251  i!=lines.end(); ++i)
1252  if (i->inhomogeneity != 0.)
1253  return true;
1254 
1255  return false;
1256 }
1257 
1258 
1259 void ConstraintMatrix::print (std::ostream &out) const
1260 {
1261  for (size_type i=0; i!=lines.size(); ++i)
1262  {
1263  // output the list of constraints as pairs of dofs and their weights
1264  if (lines[i].entries.size() > 0)
1265  {
1266  for (size_type j=0; j<lines[i].entries.size(); ++j)
1267  out << " " << lines[i].index
1268  << " " << lines[i].entries[j].first
1269  << ": " << lines[i].entries[j].second << "\n";
1270 
1271  // print out inhomogeneity.
1272  if (lines[i].inhomogeneity != 0)
1273  out << " " << lines[i].index
1274  << ": " << lines[i].inhomogeneity << "\n";
1275  }
1276  else
1277  // but also output something if the constraint simply reads
1278  // x[13]=0, i.e. where the right hand side is not a linear
1279  // combination of other dofs
1280  {
1281  if (lines[i].inhomogeneity != 0)
1282  out << " " << lines[i].index
1283  << " = " << lines[i].inhomogeneity
1284  << "\n";
1285  else
1286  out << " " << lines[i].index << " = 0\n";
1287  }
1288  }
1289 
1290  AssertThrow (out, ExcIO());
1291 }
1292 
1293 
1294 
1295 void
1296 ConstraintMatrix::write_dot (std::ostream &out) const
1297 {
1298  out << "digraph constraints {"
1299  << std::endl;
1300  for (size_type i=0; i!=lines.size(); ++i)
1301  {
1302  // same concept as in the previous function
1303  if (lines[i].entries.size() > 0)
1304  for (size_type j=0; j<lines[i].entries.size(); ++j)
1305  out << " " << lines[i].index << "->" << lines[i].entries[j].first
1306  << "; // weight: "
1307  << lines[i].entries[j].second
1308  << "\n";
1309  else
1310  out << " " << lines[i].index << "\n";
1311  }
1312  out << "}" << std::endl;
1313 }
1314 
1315 
1316 
1317 std::size_t
1319 {
1324 }
1325 
1326 
1327 
1328 void
1329 ConstraintMatrix::resolve_indices (std::vector<types::global_dof_index> &indices) const
1330 {
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)
1334  {
1335  line_ptr = get_constraint_entries(indices[i]);
1336  // if the index is constraint, the constraints indices are added to the
1337  // indices vector
1338  if (line_ptr!=nullptr)
1339  {
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);
1343  }
1344  }
1345 
1346  // keep only the unique elements
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());
1351 }
1352 
1353 
1354 
1355 // explicit instantiations
1356 //
1357 // define a list of functions for vectors and matrices, respectively, where
1358 // the vector/matrix can be replaced using a preprocessor variable
1359 // VectorType/MatrixType. note that we need a space between "VectorType" and
1360 // ">" to disambiguate ">>" when VectorType trails in an angle bracket
1361 
1362 // TODO: The way we define all the instantiations is probably not the very
1363 // best one. Try to find a better description.
1364 
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> &, \
1372  VectorType &, \
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> &, \
1378  VectorType &, \
1379  const FullMatrix<VectorType::value_type> &, \
1380  bool) const
1381 
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> &, \
1386  VectorType &, \
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> &, \
1392  VectorType &, \
1393  const FullMatrix<VectorType::value_type> &, \
1394  bool) const
1395 
1396 #ifdef DEAL_II_WITH_PETSC
1397 VECTOR_FUNCTIONS(PETScWrappers::MPI::Vector);
1398 VECTOR_FUNCTIONS(PETScWrappers::MPI::BlockVector);
1399 #endif
1400 
1401 #ifdef DEAL_II_WITH_TRILINOS
1402 PARALLEL_VECTOR_FUNCTIONS(TrilinosWrappers::MPI::Vector);
1403 PARALLEL_VECTOR_FUNCTIONS(TrilinosWrappers::MPI::BlockVector);
1404 #endif
1405 
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> &, \
1411  MatrixType &, \
1412  VectorType &, \
1413  bool , \
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> &, \
1420  MatrixType &, \
1421  Vector<VectorScalar> &, \
1422  bool , \
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> &, \
1429  MatrixType &, \
1430  VectorType &, \
1431  bool , \
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> &, \
1438  MatrixType &, \
1439  Vector<MatrixType::value_type> &, \
1440  bool , \
1441  std::integral_constant<bool, true>) const
1442 
1443 MATRIX_FUNCTIONS(FullMatrix<double>,double);
1444 MATRIX_FUNCTIONS(FullMatrix<float>,float);
1445 MATRIX_FUNCTIONS(FullMatrix<double>,std::complex<double>);
1446 MATRIX_FUNCTIONS(FullMatrix<std::complex<double> >,std::complex<double>);
1447 
1448 MATRIX_FUNCTIONS(SparseMatrix<double>,double);
1449 MATRIX_FUNCTIONS(SparseMatrix<float>,float);
1450 MATRIX_FUNCTIONS(SparseMatrix<double>,std::complex<double>);
1451 MATRIX_FUNCTIONS(SparseMatrix<float>,std::complex<float>);
1452 MATRIX_FUNCTIONS(SparseMatrix<std::complex<double> >,std::complex<double>);
1453 MATRIX_FUNCTIONS(SparseMatrix<std::complex<float> >,std::complex<float>);
1454 
1455 MATRIX_FUNCTIONS(SparseMatrixEZ<double>,double);
1456 MATRIX_FUNCTIONS(SparseMatrixEZ<float>,float);
1457 MATRIX_FUNCTIONS(ChunkSparseMatrix<double>,double);
1458 MATRIX_FUNCTIONS(ChunkSparseMatrix<float>,float);
1459 
1460 
1461 BLOCK_MATRIX_FUNCTIONS(BlockSparseMatrix<double>);
1462 BLOCK_MATRIX_FUNCTIONS(BlockSparseMatrix<float>);
1463 BLOCK_MATRIX_VECTOR_FUNCTIONS(BlockSparseMatrix<double>, BlockVector<double>);
1464 BLOCK_MATRIX_VECTOR_FUNCTIONS(BlockSparseMatrix<float>, BlockVector<float>);
1465 
1466 // BLOCK_MATRIX_FUNCTIONS(BlockSparseMatrixEZ<double>);
1467 // BLOCK_MATRIX_VECTOR_FUNCTIONS(BlockSparseMatrixEZ<float>, Vector<float>);
1468 
1469 #ifdef DEAL_II_WITH_PETSC
1470 MATRIX_FUNCTIONS(PETScWrappers::SparseMatrix,PetscScalar);
1471 MATRIX_FUNCTIONS(PETScWrappers::MPI::SparseMatrix,PetscScalar);
1472 BLOCK_MATRIX_FUNCTIONS(PETScWrappers::MPI::BlockSparseMatrix);
1476 #endif
1477 
1478 #ifdef DEAL_II_WITH_TRILINOS
1479 MATRIX_FUNCTIONS(TrilinosWrappers::SparseMatrix,double);
1480 BLOCK_MATRIX_FUNCTIONS(TrilinosWrappers::BlockSparseMatrix);
1483 #endif
1484 
1485 
1486 #define SPARSITY_FUNCTIONS(SparsityPatternType) \
1487  template void ConstraintMatrix::add_entries_local_to_global<SparsityPatternType> ( \
1488  const std::vector<ConstraintMatrix::size_type> &, \
1489  SparsityPatternType &, \
1490  const bool, \
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 &, \
1497  const bool, \
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 &, \
1503  const bool, \
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 &, \
1510  const bool, \
1511  const Table<2,bool> &) const
1512 
1513 SPARSITY_FUNCTIONS(SparsityPattern);
1514 SPARSITY_FUNCTIONS(DynamicSparsityPattern);
1515 BLOCK_SPARSITY_FUNCTIONS(BlockSparsityPattern);
1516 BLOCK_SPARSITY_FUNCTIONS(BlockDynamicSparsityPattern);
1517 
1518 #ifdef DEAL_II_WITH_TRILINOS
1519 SPARSITY_FUNCTIONS(TrilinosWrappers::SparsityPattern);
1520 BLOCK_SPARSITY_FUNCTIONS(TrilinosWrappers::BlockSparsityPattern);
1521 #endif
1522 
1523 
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> &, \
1535  MatrixType &) const
1536 
1537 ONLY_MATRIX_FUNCTIONS(FullMatrix<float>);
1538 ONLY_MATRIX_FUNCTIONS(FullMatrix<double>);
1539 ONLY_MATRIX_FUNCTIONS(SparseMatrix<float>);
1540 ONLY_MATRIX_FUNCTIONS(SparseMatrix<double>);
1541 ONLY_MATRIX_FUNCTIONS(MatrixBlock<SparseMatrix<float> >);
1542 ONLY_MATRIX_FUNCTIONS(MatrixBlock<SparseMatrix<double> >);
1543 ONLY_MATRIX_FUNCTIONS(BlockSparseMatrix<float>);
1544 ONLY_MATRIX_FUNCTIONS(BlockSparseMatrix<double>);
1545 
1546 #ifdef DEAL_II_WITH_TRILINOS
1547 ONLY_MATRIX_FUNCTIONS(TrilinosWrappers::SparseMatrix);
1548 ONLY_MATRIX_FUNCTIONS(TrilinosWrappers::BlockSparseMatrix);
1549 #endif
1550 
1551 #ifdef DEAL_II_WITH_PETSC
1552 ONLY_MATRIX_FUNCTIONS(PETScWrappers::SparseMatrix);
1553 ONLY_MATRIX_FUNCTIONS(PETScWrappers::MPI::SparseMatrix);
1554 ONLY_MATRIX_FUNCTIONS(PETScWrappers::MPI::BlockSparseMatrix);
1555 #endif
1556 
1557 #include "constraint_matrix.inst"
1558 
1559 // allocate scratch data. Cannot use the generic template instantiation
1560 // because we need to provide an initializer object of type
1561 // internals::ConstraintMatrixData<Number> that can be passed to the
1562 // constructor of scratch_data (it won't allow one to be constructed in place).
1563 namespace internals
1564 {
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)
1569 
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
1577 }
1578 
1579 
1580 DEAL_II_NAMESPACE_CLOSE
const types::global_dof_index invalid_size_type
Definition: types.h:182
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)
boost::iterator_range< const_iterator > LineRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
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)
Definition: index_set.h:1619
iterator end() const
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
types::global_dof_index size_type
std::vector< std::pair< size_type, double > > Entries
static ::ExceptionBase & ExcLineInexistant(size_type arg1)
size_type n_cols() const
bool is_compressed() const
SparsityPatternType & block(const size_type row, const size_type column)
size_type size() const
Definition: index_set.h:1575
size_type n_rows() 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)
Definition: index_set.cc:288
#define Assert(cond, exc)
Definition: exceptions.h:1142
void print(std::ostream &out) const
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1811
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)
Definition: vector.h:1312
void distribute(VectorType &vec) const
void compress() const
Definition: index_set.h:1584
void write_dot(std::ostream &) const
void copy_from(const ConstraintMatrix &other)
std::atomic< unsigned int > counter
Definition: subscriptor.h:203
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
Definition: index_set.h:1645
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)
iterator begin() const
const LineRange get_lines() const
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) 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()