Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
index_set.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2005 - 2020 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.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
18#include <deal.II/base/mpi.h>
19
20#include <vector>
21
22#ifdef DEAL_II_WITH_TRILINOS
23# ifdef DEAL_II_WITH_MPI
24# include <Epetra_MpiComm.h>
25# endif
26# include <Epetra_Map.h>
27# include <Epetra_SerialComm.h>
28#endif
29
31
32
33
34#ifdef DEAL_II_WITH_TRILINOS
35
36// the 64-bit path uses a few different names, so put that into a separate
37// implementation
38
39# ifdef DEAL_II_WITH_64BIT_INDICES
40
41IndexSet::IndexSet(const Epetra_BlockMap &map)
42 : is_compressed(true)
43 , index_space_size(1 + map.MaxAllGID64())
44 , largest_range(numbers::invalid_unsigned_int)
45{
46 Assert(map.MinAllGID64() == 0,
48 "The Epetra_BlockMap does not contain the global index 0, "
49 "which means some entries are not present on any processor."));
50
51 // For a contiguous map, we do not need to go through the whole data...
52 if (map.LinearMap())
53 add_range(size_type(map.MinMyGID64()), size_type(map.MaxMyGID64() + 1));
54 else
55 {
56 const size_type n_indices = map.NumMyElements();
57 size_type * indices =
58 reinterpret_cast<size_type *>(map.MyGlobalElements64());
59 add_indices(indices, indices + n_indices);
60 }
61 compress();
62}
63
64# else
65
66// this is the standard 32-bit implementation
67
68IndexSet::IndexSet(const Epetra_BlockMap &map)
69 : is_compressed(true)
70 , index_space_size(1 + map.MaxAllGID())
71 , largest_range(numbers::invalid_unsigned_int)
72{
73 Assert(map.MinAllGID() == 0,
75 "The Epetra_BlockMap does not contain the global index 0, "
76 "which means some entries are not present on any processor."));
77
78 // For a contiguous map, we do not need to go through the whole data...
79 if (map.LinearMap())
80 add_range(size_type(map.MinMyGID()), size_type(map.MaxMyGID() + 1));
81 else
82 {
83 const size_type n_indices = map.NumMyElements();
84 unsigned int * indices =
85 reinterpret_cast<unsigned int *>(map.MyGlobalElements());
86 add_indices(indices, indices + n_indices);
87 }
88 compress();
89}
90
91# endif
92
93#endif // ifdef DEAL_II_WITH_TRILINOS
94
95
96
97void
99{
100 // we will, in the following, modify mutable variables. this can only
101 // work in multithreaded applications if we lock the data structures
102 // via a mutex, so that users can call 'const' functions from threads
103 // in parallel (and these 'const' functions can then call compress()
104 // which itself calls the current function)
105 std::lock_guard<std::mutex> lock(compress_mutex);
106
107 // see if any of the contiguous ranges can be merged. do not use
108 // std::vector::erase in-place as it is quadratic in the number of
109 // ranges. since the ranges are sorted by their first index, determining
110 // overlap isn't all that hard
111 std::vector<Range>::iterator store = ranges.begin();
112 for (std::vector<Range>::iterator i = ranges.begin(); i != ranges.end();)
113 {
114 std::vector<Range>::iterator next = i;
115 ++next;
116
117 size_type first_index = i->begin;
118 size_type last_index = i->end;
119
120 // see if we can merge any of the following ranges
121 while (next != ranges.end() && (next->begin <= last_index))
122 {
123 last_index = std::max(last_index, next->end);
124 ++next;
125 }
126 i = next;
127
128 // store the new range in the slot we last occupied
129 *store = Range(first_index, last_index);
130 ++store;
131 }
132 // use a compact array with exactly the right amount of storage
133 if (store != ranges.end())
134 {
135 std::vector<Range> new_ranges(ranges.begin(), store);
136 ranges.swap(new_ranges);
137 }
138
139 // now compute indices within set and the range with most elements
140 size_type next_index = 0, largest_range_size = 0;
141 for (std::vector<Range>::iterator i = ranges.begin(); i != ranges.end(); ++i)
142 {
143 Assert(i->begin < i->end, ExcInternalError());
144
145 i->nth_index_in_set = next_index;
146 next_index += (i->end - i->begin);
147 if (i->end - i->begin > largest_range_size)
148 {
149 largest_range_size = i->end - i->begin;
150 largest_range = i - ranges.begin();
151 }
152 }
153 is_compressed = true;
154
155 // check that next_index is correct. needs to be after the previous
156 // statement because we otherwise will get into an endless loop
157 Assert(next_index == n_elements(), ExcInternalError());
158}
159
160
161
162#ifndef DOXYGEN
164{
165 Assert(size() == is.size(), ExcDimensionMismatch(size(), is.size()));
166
167 compress();
168 is.compress();
169
170 std::vector<Range>::const_iterator r1 = ranges.begin(),
171 r2 = is.ranges.begin();
172 IndexSet result(size());
173
174 while ((r1 != ranges.end()) && (r2 != is.ranges.end()))
175 {
176 // if r1 and r2 do not overlap at all, then move the pointer that sits
177 // to the left of the other up by one
178 if (r1->end <= r2->begin)
179 ++r1;
180 else if (r2->end <= r1->begin)
181 ++r2;
182 else
183 {
184 // the ranges must overlap somehow
185 Assert(((r1->begin <= r2->begin) && (r1->end > r2->begin)) ||
186 ((r2->begin <= r1->begin) && (r2->end > r1->begin)),
188
189 // add the overlapping range to the result
190 result.add_range(std::max(r1->begin, r2->begin),
191 std::min(r1->end, r2->end));
192
193 // now move that iterator that ends earlier one up. note that it has
194 // to be this one because a subsequent range may still have a chance
195 // of overlapping with the range that ends later
196 if (r1->end <= r2->end)
197 ++r1;
198 else
199 ++r2;
200 }
201 }
202
203 result.compress();
204 return result;
205}
206#endif
207
208
209
212{
213 Assert(begin <= end,
214 ExcMessage("End index needs to be larger or equal to begin index!"));
215 Assert(end <= size(), ExcMessage("Given range exceeds index set dimension"));
216
217 IndexSet result(end - begin);
218 std::vector<Range>::const_iterator r1 = ranges.begin();
219
220 while (r1 != ranges.end())
221 {
222 if ((r1->end > begin) && (r1->begin < end))
223 {
224 result.add_range(std::max(r1->begin, begin) - begin,
225 std::min(r1->end, end) - begin);
226 }
227 else if (r1->begin >= end)
228 break;
229
230 ++r1;
231 }
232
233 result.compress();
234 return result;
235}
236
237std::vector<IndexSet>
239 const std::vector<types::global_dof_index> &n_indices_per_block) const
240{
241 std::vector<IndexSet> partitioned;
242 const unsigned int n_blocks = n_indices_per_block.size();
243
244 partitioned.reserve(n_blocks);
245 types::global_dof_index start = 0;
247 for (const auto n_block_indices : n_indices_per_block)
248 {
249 partitioned.push_back(this->get_view(start, start + n_block_indices));
250 start += n_block_indices;
251 sum += partitioned.back().size();
252 }
253 AssertDimension(sum, this->size());
254 return partitioned;
255}
256
257void
259{
260 compress();
261 other.compress();
262 is_compressed = false;
263
264
265 // we save new ranges to be added to our IndexSet in an temporary vector and
266 // add all of them in one go at the end.
267 std::vector<Range> new_ranges;
268
269 std::vector<Range>::iterator own_it = ranges.begin();
270 std::vector<Range>::iterator other_it = other.ranges.begin();
271
272 while (own_it != ranges.end() && other_it != other.ranges.end())
273 {
274 // advance own iterator until we get an overlap
275 if (own_it->end <= other_it->begin)
276 {
277 ++own_it;
278 continue;
279 }
280 // we are done with other_it, so advance
281 if (own_it->begin >= other_it->end)
282 {
283 ++other_it;
284 continue;
285 }
286
287 // Now own_it and other_it overlap. First save the part of own_it that
288 // is before other_it (if not empty).
289 if (own_it->begin < other_it->begin)
290 {
291 Range r(own_it->begin, other_it->begin);
292 r.nth_index_in_set = 0; // fix warning of unused variable
293 new_ranges.push_back(r);
294 }
295 // change own_it to the sub range behind other_it. Do not delete own_it
296 // in any case. As removal would invalidate iterators, we just shrink
297 // the range to an empty one.
298 own_it->begin = other_it->end;
299 if (own_it->begin > own_it->end)
300 {
301 own_it->begin = own_it->end;
302 ++own_it;
303 }
304
305 // continue without advancing iterators, the right one will be advanced
306 // next.
307 }
308
309 // Now delete all empty ranges we might
310 // have created.
311 for (std::vector<Range>::iterator it = ranges.begin(); it != ranges.end();)
312 {
313 if (it->begin >= it->end)
314 it = ranges.erase(it);
315 else
316 ++it;
317 }
318
319 // done, now add the temporary ranges
320 const std::vector<Range>::iterator end = new_ranges.end();
321 for (std::vector<Range>::iterator it = new_ranges.begin(); it != end; ++it)
322 add_range(it->begin, it->end);
323
324 compress();
325}
326
327
328
331{
332 IndexSet set(this->size() * other.size());
333 for (const auto el : *this)
334 set.add_indices(other, el * other.size());
335 set.compress();
336 return set;
337}
338
339
340
343{
344 Assert(is_empty() == false,
346 "pop_back() failed, because this IndexSet contains no entries."));
347
348 const size_type index = ranges.back().end - 1;
349 --ranges.back().end;
350
351 if (ranges.back().begin == ranges.back().end)
352 ranges.pop_back();
353
354 return index;
355}
356
357
358
361{
362 Assert(is_empty() == false,
364 "pop_front() failed, because this IndexSet contains no entries."));
365
366 const size_type index = ranges.front().begin;
367 ++ranges.front().begin;
368
369 if (ranges.front().begin == ranges.front().end)
370 ranges.erase(ranges.begin());
371
372 // We have to set this in any case, because nth_index_in_set is no longer
373 // up to date for all but the first range
374 is_compressed = false;
375
376 return index;
377}
378
379
380
381void
382IndexSet::add_indices(const IndexSet &other, const size_type offset)
383{
384 if ((this == &other) && (offset == 0))
385 return;
386
387 if (other.ranges.size() != 0)
388 {
389 AssertIndexRange(other.ranges.back().end - 1, index_space_size);
390 }
391
392 compress();
393 other.compress();
394
395 std::vector<Range>::const_iterator r1 = ranges.begin(),
396 r2 = other.ranges.begin();
397
398 std::vector<Range> new_ranges;
399 // just get the start and end of the ranges right in this method, everything
400 // else will be done in compress()
401 while (r1 != ranges.end() || r2 != other.ranges.end())
402 {
403 // the two ranges do not overlap or we are at the end of one of the
404 // ranges
405 if (r2 == other.ranges.end() ||
406 (r1 != ranges.end() && r1->end < (r2->begin + offset)))
407 {
408 new_ranges.push_back(*r1);
409 ++r1;
410 }
411 else if (r1 == ranges.end() || (r2->end + offset) < r1->begin)
412 {
413 new_ranges.emplace_back(r2->begin + offset, r2->end + offset);
414 ++r2;
415 }
416 else
417 {
418 // ok, we do overlap, so just take the combination of the current
419 // range (do not bother to merge with subsequent ranges)
420 Range next(std::min(r1->begin, r2->begin + offset),
421 std::max(r1->end, r2->end + offset));
422 new_ranges.push_back(next);
423 ++r1;
424 ++r2;
425 }
426 }
427 ranges.swap(new_ranges);
428
429 is_compressed = false;
430 compress();
431}
432
433
434
435void
436IndexSet::write(std::ostream &out) const
437{
438 compress();
439 out << size() << " ";
440 out << ranges.size() << std::endl;
441 std::vector<Range>::const_iterator r = ranges.begin();
442 for (; r != ranges.end(); ++r)
443 {
444 out << r->begin << " " << r->end << std::endl;
445 }
446}
447
448
449
450void
451IndexSet::read(std::istream &in)
452{
453 AssertThrow(in, ExcIO());
454
455 size_type s;
456 unsigned int n_ranges;
457
458 in >> s >> n_ranges;
459 ranges.clear();
460 set_size(s);
461 for (unsigned int i = 0; i < n_ranges; ++i)
462 {
463 AssertThrow(in, ExcIO());
464
465 size_type b, e;
466 in >> b >> e;
467 add_range(b, e);
468 }
469}
470
471
472void
473IndexSet::block_write(std::ostream &out) const
474{
475 AssertThrow(out, ExcIO());
476 out.write(reinterpret_cast<const char *>(&index_space_size),
477 sizeof(index_space_size));
478 std::size_t n_ranges = ranges.size();
479 out.write(reinterpret_cast<const char *>(&n_ranges), sizeof(n_ranges));
480 if (ranges.empty() == false)
481 out.write(reinterpret_cast<const char *>(&*ranges.begin()),
482 ranges.size() * sizeof(Range));
483 AssertThrow(out, ExcIO());
484}
485
486void
487IndexSet::block_read(std::istream &in)
488{
490 std::size_t n_ranges;
491 in.read(reinterpret_cast<char *>(&size), sizeof(size));
492 in.read(reinterpret_cast<char *>(&n_ranges), sizeof(n_ranges));
493 // we have to clear ranges first
494 ranges.clear();
495 set_size(size);
496 ranges.resize(n_ranges, Range(0, 0));
497 if (n_ranges)
498 in.read(reinterpret_cast<char *>(&*ranges.begin()),
499 ranges.size() * sizeof(Range));
500
501 do_compress(); // needed so that largest_range can be recomputed
502}
503
504
505
506void
507IndexSet::fill_index_vector(std::vector<size_type> &indices) const
508{
509 compress();
510
511 indices.clear();
512 indices.reserve(n_elements());
513
514 for (const auto &range : ranges)
515 for (size_type entry = range.begin; entry < range.end; ++entry)
516 indices.push_back(entry);
517
518 Assert(indices.size() == n_elements(), ExcInternalError());
519}
520
521
522
523#ifdef DEAL_II_WITH_TRILINOS
524# ifdef DEAL_II_TRILINOS_WITH_TPETRA
525
526Tpetra::Map<int, types::global_dof_index>
528 const bool overlapping) const
529{
530 compress();
531 (void)communicator;
532
533# ifdef DEBUG
534 if (!overlapping)
535 {
537 Utilities::MPI::sum(n_elements(), communicator);
539 ExcMessage("You are trying to create an Tpetra::Map object "
540 "that partitions elements of an index set "
541 "between processors. However, the union of the "
542 "index sets on different processors does not "
543 "contain all indices exactly once: the sum of "
544 "the number of entries the various processors "
545 "want to store locally is " +
547 " whereas the total size of the object to be "
548 "allocated is " +
550 ". In other words, there are "
551 "either indices that are not spoken for "
552 "by any processor, or there are indices that are "
553 "claimed by multiple processors."));
554 }
555# endif
556
557 // Find out if the IndexSet is ascending and 1:1. This corresponds to a
558 // linear Tpetra::Map. Overlapping IndexSets are never 1:1.
559 const bool linear =
560 overlapping ? false : is_ascending_and_one_to_one(communicator);
561 if (linear)
562 return Tpetra::Map<int, types::global_dof_index>(
563 size(),
564 n_elements(),
565 0,
566# ifdef DEAL_II_WITH_MPI
567 Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
568# else
569 Teuchos::rcp(new Teuchos::Comm<int>())
570# endif
571 );
572 else
573 {
574 std::vector<size_type> indices;
575 fill_index_vector(indices);
576 std::vector<types::global_dof_index> int_indices(indices.size());
577 std::copy(indices.begin(), indices.end(), int_indices.begin());
578 const Teuchos::ArrayView<types::global_dof_index> arr_view(int_indices);
579 return Tpetra::Map<int, types::global_dof_index>(
580 size(),
581 arr_view,
582 0,
583# ifdef DEAL_II_WITH_MPI
584 Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
585# else
586 Teuchos::rcp(new Teuchos::Comm<int>())
587# endif
588 );
589 }
590}
591# endif
592
593
594
595Epetra_Map
597 const bool overlapping) const
598{
599 compress();
600 (void)communicator;
601
602# ifdef DEBUG
603 if (!overlapping)
604 {
606 Utilities::MPI::sum(n_elements(), communicator);
608 ExcMessage("You are trying to create an Epetra_Map object "
609 "that partitions elements of an index set "
610 "between processors. However, the union of the "
611 "index sets on different processors does not "
612 "contain all indices exactly once: the sum of "
613 "the number of entries the various processors "
614 "want to store locally is " +
616 " whereas the total size of the object to be "
617 "allocated is " +
619 ". In other words, there are "
620 "either indices that are not spoken for "
621 "by any processor, or there are indices that are "
622 "claimed by multiple processors."));
623 }
624# endif
625
626 // Find out if the IndexSet is ascending and 1:1. This corresponds to a
627 // linear EpetraMap. Overlapping IndexSets are never 1:1.
628 const bool linear =
629 overlapping ? false : is_ascending_and_one_to_one(communicator);
630
631 if (linear)
632 return Epetra_Map(TrilinosWrappers::types::int_type(size()),
634 0,
635# ifdef DEAL_II_WITH_MPI
636 Epetra_MpiComm(communicator)
637# else
638 Epetra_SerialComm()
639# endif
640 );
641 else
642 {
643 std::vector<size_type> indices;
644 fill_index_vector(indices);
645 return Epetra_Map(
648 (n_elements() > 0 ?
649 reinterpret_cast<TrilinosWrappers::types::int_type *>(
650 indices.data()) :
651 nullptr),
652 0,
653# ifdef DEAL_II_WITH_MPI
654 Epetra_MpiComm(communicator)
655# else
656 Epetra_SerialComm()
657# endif
658 );
659 }
660}
661#endif
662
663
664
665bool
667{
668 // If the sum of local elements does not add up to the total size,
669 // the IndexSet can't be complete.
671 Utilities::MPI::sum(n_elements(), communicator);
672 if (n_global_elements != size())
673 return false;
674
675 if (n_global_elements == 0)
676 return true;
677
678#ifdef DEAL_II_WITH_MPI
679 // Non-contiguous IndexSets can't be linear.
680 const bool all_contiguous =
681 (Utilities::MPI::min(is_contiguous() ? 1 : 0, communicator) == 1);
682 if (!all_contiguous)
683 return false;
684
685 bool is_globally_ascending = true;
686 // we know that there is only one interval
687 types::global_dof_index first_local_dof = (n_elements() > 0) ?
688 *(begin_intervals()->begin()) :
690
691 const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
692 const std::vector<types::global_dof_index> global_dofs =
693 Utilities::MPI::gather(communicator, first_local_dof, 0);
694
695 if (my_rank == 0)
696 {
697 // find out if the received std::vector is ascending
698 types::global_dof_index index = 0;
699 while (global_dofs[index] == numbers::invalid_dof_index)
700 ++index;
701 types::global_dof_index old_dof = global_dofs[index++];
702 for (; index < global_dofs.size(); ++index)
703 {
704 const types::global_dof_index new_dof = global_dofs[index];
705 if (new_dof != numbers::invalid_dof_index)
706 {
707 if (new_dof <= old_dof)
708 {
709 is_globally_ascending = false;
710 break;
711 }
712 else
713 old_dof = new_dof;
714 }
715 }
716 }
717
718 // now broadcast the result
719 int is_ascending = is_globally_ascending ? 1 : 0;
720 int ierr = MPI_Bcast(&is_ascending, 1, MPI_INT, 0, communicator);
721 AssertThrowMPI(ierr);
722
723 return (is_ascending == 1);
724#else
725 return true;
726#endif // DEAL_II_WITH_MPI
727}
728
729
730
731std::size_t
733{
737 sizeof(compress_mutex));
738}
739
740
741
ElementIterator begin() const
Definition: index_set.h:1074
size_type largest_range
Definition: index_set.h:979
bool is_contiguous() const
Definition: index_set.h:1815
void do_compress() const
Definition: index_set.cc:98
size_type size() const
Definition: index_set.h:1634
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:507
bool is_empty() const
Definition: index_set.h:1824
std::vector< IndexSet > split_by_block(const std::vector< types::global_dof_index > &n_indices_per_block) const
Definition: index_set.cc:238
size_type n_elements() const
Definition: index_set.h:1832
ElementIterator begin() const
Definition: index_set.h:1518
size_type pop_front()
Definition: index_set.cc:360
void set_size(const size_type size)
Definition: index_set.h:1622
void read(std::istream &in)
Definition: index_set.cc:451
Tpetra::Map< int, types::global_dof_index > make_tpetra_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:527
IndexSet tensor_product(const IndexSet &other) const
Definition: index_set.cc:330
void write(std::ostream &out) const
Definition: index_set.cc:436
void block_read(std::istream &in)
Definition: index_set.cc:487
bool is_compressed
Definition: index_set.h:962
IntervalIterator begin_intervals() const
Definition: index_set.h:1589
std::vector< Range > ranges
Definition: index_set.h:952
void subtract_set(const IndexSet &other)
Definition: index_set.cc:258
ElementIterator end() const
Definition: index_set.h:1580
Threads::Mutex compress_mutex
Definition: index_set.h:985
size_type index_space_size
Definition: index_set.h:968
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:596
void block_write(std::ostream &out) const
Definition: index_set.cc:473
IndexSet get_view(const size_type begin, const size_type end) const
Definition: index_set.cc:211
void add_range(const size_type begin, const size_type end)
Definition: index_set.h:1673
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:666
std::size_t memory_consumption() const
Definition: index_set.cc:732
void compress() const
Definition: index_set.h:1642
size_type pop_back()
Definition: index_set.cc:342
types::global_dof_index size_type
Definition: index_set.h:83
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1703
IndexSet operator&(const IndexSet &is) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
__global__ void set(Number *val, const Number s, const size_type N)
static ::ExceptionBase & ExcIO()
#define Assert(cond, exc)
Definition: exceptions.h:1465
std::string to_string(const T &t)
Definition: patterns.h:2329
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1746
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
types::global_dof_index size_type
Definition: cuda_kernels.h:45
std::enable_if< IsBlockVector< VectorType >::value, unsignedint >::type n_blocks(const VectorType &vector)
Definition: operators.h:49
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
T min(const T &t, const MPI_Comm &mpi_communicator)
std::vector< T > gather(const MPI_Comm &comm, const T &object_to_send, const unsigned int root_process=0)
T sum(const T &t, const MPI_Comm &mpi_communicator)
std::string compress(const std::string &input)
Definition: utilities.cc:392
void copy(const T *begin, const T *end, U *dest)
static const unsigned int invalid_unsigned_int
Definition: types.h:196
const types::global_dof_index invalid_dof_index
Definition: types.h:211
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
size_type nth_index_in_set
Definition: index_set.h:880