opm-common
Loading...
Searching...
No Matches
CSRGraphFromCoordinates_impl.hpp
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 Statoil ASA.
4 Copyright 2022 Equinor ASA
5
6 This file is part of the Open Porous Media Project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#include <algorithm>
23#include <cassert>
24#include <exception>
25#include <iterator>
26#include <optional>
27#include <set>
28#include <stdexcept>
29#include <string>
30#include <unordered_map>
31#include <utility>
32#include <vector>
33
34// ---------------------------------------------------------------------
35// Class Opm::utility::CSRGraphFromCoordinates::Connections
36// ---------------------------------------------------------------------
37
38template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
39void
41Connections::add(const VertexID v1, const VertexID v2)
42{
43 this->i_.push_back(v1);
44 this->j_.push_back(v2);
45
46 this->max_i_ = std::max(this->max_i_.value_or(BaseVertexID{}), this->i_.back());
47 this->max_j_ = std::max(this->max_j_.value_or(BaseVertexID{}), this->j_.back());
48}
49
50template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
51void
53Connections::add(VertexID maxRowIdx,
54 VertexID maxColIdx,
55 const Neighbours& rows,
56 const Neighbours& cols)
57{
58 if (cols.size() != rows.size()) {
59 throw std::invalid_argument {
60 "Coordinate format column index table size does not match "
61 "row index table size"
62 };
63 }
64
65 this->i_.insert(this->i_.end(), rows .begin(), rows .end());
66 this->j_.insert(this->j_.end(), cols .begin(), cols .end());
67
68 this->max_i_ = std::max(this->max_i_.value_or(BaseVertexID{}), maxRowIdx);
69 this->max_j_ = std::max(this->max_j_.value_or(BaseVertexID{}), maxColIdx);
70}
71
72template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
73void
76{
77 this->j_.clear();
78 this->i_.clear();
79
80 this->max_i_.reset();
81 this->max_j_.reset();
82}
83
84template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
85bool
88{
89 return this->i_.empty();
90}
91
92template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
93bool
96{
97 return this->i_.size() == this->j_.size();
98}
99
100template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
101std::optional<typename Opm::utility::CSRGraphFromCoordinates<VertexID, TrackCompressedIdx, PermitSelfConnections>::BaseVertexID>
104{
105 return this->max_i_;
106}
107
108template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
109std::optional<typename Opm::utility::CSRGraphFromCoordinates<VertexID, TrackCompressedIdx, PermitSelfConnections>::BaseVertexID>
112{
113 return this->max_j_;
114}
115
116template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
120{
121 return this->i_.size();
122}
123
124template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
128{
129 return this->i_;
130}
131
132template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
136{
137 return this->j_;
138}
139
140template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
141VertexID
143Connections::findMergedVertexID(VertexID v, const std::unordered_map<VertexID, VertexID>& vertex_merges) const
144{
145 // If found in the map, return the merged vertex ID
146 // Otherwise, return the original vertex ID unchanged
147 // vertex_merges is fully resolved from our disjoint set union structure
148 auto it = vertex_merges.find(v);
149 if (it != vertex_merges.end()) {
150 return it->second;
151 } else {
152 return v;
153 }
154}
155
156
157template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
158std::unordered_map<VertexID, VertexID>
160Connections::applyVertexMerges(const std::unordered_map<VertexID, VertexID>& vertex_merges)
161{
162
163 // Find the maximum vertex ID in the original connections
164 VertexID max_original_vertex_id = std::max(max_i_.value_or(BaseVertexID{}), max_j_.value_or(BaseVertexID{}));
165
166 // Apply merges to all connections directly in one pass
167 // We can leverage that vertex_merges is already fully resolved from our disjoint set union
168 // structure, so we don't need nested lookups
169 std::ranges::transform(i_, i_.begin(),
170 [this, &vertex_merges](VertexID v) { return findMergedVertexID(v, vertex_merges); });
171
172 std::ranges::transform(j_, j_.begin(),
173 [this, &vertex_merges](VertexID v) { return findMergedVertexID(v, vertex_merges); });
174
175 // Remove self-connections if required using erase-remove idiom
176 if constexpr (!PermitSelfConnections) {
177 auto write_pos = 0*i_.size();
178 for (auto read_pos = 0*i_.size(); read_pos < i_.size(); ++read_pos) {
179 if (i_[read_pos] != j_[read_pos]) {
180 if (write_pos != read_pos) {
181 i_[write_pos] = i_[read_pos];
182 j_[write_pos] = j_[read_pos];
183 }
184 ++write_pos;
185 }
186 }
187 i_.resize(write_pos);
188 j_.resize(write_pos);
189 }
190
191 // Create compact vertex numbering
192 std::set<VertexID> sorted_unique_vertices;
193 sorted_unique_vertices.insert(i_.begin(), i_.end());
194 sorted_unique_vertices.insert(j_.begin(), j_.end());
195
196 // Generate direct mapping from old to compact vertex IDs
197 std::unordered_map<VertexID, VertexID> vertex_map;
198 vertex_map.reserve(sorted_unique_vertices.size());
199 VertexID new_id = 0;
200 for (auto& v : sorted_unique_vertices) {
201 vertex_map.emplace(v, new_id++);
202 }
203
204 // Update max indices
205 this->max_i_ = this->max_j_ = new_id - 1;
206
207 // Remap all vertices to compact IDs
208 auto remap = [&vertex_map](auto v) {
209 // Using at() is appropriate here since we know all values from i_ and j_ are in vertex_map
210 return vertex_map.at(v);
211 };
212
213 std::ranges::transform(i_, i_.begin(), remap);
214 std::ranges::transform(j_, j_.begin(), remap);
215
216 // Build final mapping (original -> compact ID)
217 std::unordered_map<VertexID, VertexID> final_mapping;
218 final_mapping.reserve(max_original_vertex_id + 1);
219
220 for (VertexID vertex = 0; vertex <= max_original_vertex_id; ++vertex) {
221 auto merged_id = findMergedVertexID(vertex, vertex_merges);
222 final_mapping.emplace(vertex, vertex_map.at(merged_id));
223 }
224 return final_mapping;
225}
226
227// =====================================================================
228
229// ---------------------------------------------------------------------
230// Class Opm::utility::CSRGraphFromCoordinates::CSR
231// ---------------------------------------------------------------------
232
233template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
234void
236CSR::merge(const Connections& conns,
237 const Offset maxNumVertices,
238 const bool expandExistingIdxMap)
239{
240 const auto maxRow = conns.maxRow();
241
242 if (maxRow.has_value() &&
243 (static_cast<Offset>(*maxRow) >= maxNumVertices))
244 {
245 throw std::invalid_argument {
246 "Number of vertices in input graph (" +
247 std::to_string(*maxRow) + ") "
248 "exceeds maximum graph size implied by explicit size of "
249 "adjacency matrix (" + std::to_string(maxNumVertices) + ')'
250 };
251 }
252
253 this->assemble(conns.rowIndices(), conns.columnIndices(),
254 maxRow.value_or(BaseVertexID{0}),
255 conns.maxCol().value_or(BaseVertexID{0}),
256 expandExistingIdxMap);
257
258 this->compress(maxNumVertices);
259}
260
261template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
264CSR::numRows() const
265{
266 return this->startPointers().empty()
267 ? 0 : this->startPointers().size() - 1;
268}
269
270template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
271typename Opm::utility::CSRGraphFromCoordinates<VertexID, TrackCompressedIdx, PermitSelfConnections>::BaseVertexID
273CSR::maxRowID() const
274{
275 return this->numRows_ - 1;
276}
277
278template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
279typename Opm::utility::CSRGraphFromCoordinates<VertexID, TrackCompressedIdx, PermitSelfConnections>::BaseVertexID
281CSR::maxColID() const
282{
283 return this->numCols_ - 1;
284}
285
286template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
289CSR::startPointers() const
290{
291 return this->ia_;
292}
293
294template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
297CSR::columnIndices() const
298{
299 return this->ja_;
300}
301
302template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
306{
307 auto rowIdx = Neighbours{};
308
309 if (this->ia_.empty()) {
310 return rowIdx;
311 }
312
313 rowIdx.reserve(this->ia_.back());
314
315 auto row = BaseVertexID{};
316
317 const auto m = this->ia_.size() - 1;
318 for (auto i = 0*m; i < m; ++i, ++row) {
319 const auto n = this->ia_[i + 1] - this->ia_[i + 0];
320
321 rowIdx.insert(rowIdx.end(), n, row);
322 }
323
324 return rowIdx;
325}
326
327template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
328void
331{
332 this->ia_.clear();
333 this->ja_.clear();
334
335 if constexpr (TrackCompressedIdx) {
336 this->compressedIdx_.clear();
337 }
338
339 this->numRows_ = 0;
340 this->numCols_ = 0;
341}
342
343template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
344void
346CSR::assemble(const Neighbours& rows,
347 const Neighbours& cols,
348 const BaseVertexID maxRowIdx,
349 const BaseVertexID maxColIdx,
350 [[maybe_unused]] const bool expandExistingIdxMap)
351{
352 [[maybe_unused]] auto compressedIdx = this->compressedIdx_;
353 [[maybe_unused]] const auto numOrigNNZ = this->ja_.size();
354
355 auto i = this->coordinateFormatRowIndices();
356 i.insert(i.end(), rows.begin(), rows.end());
357
358 auto j = this->ja_;
359 j.insert(j.end(), cols.begin(), cols.end());
360
361 // Use the number of unique vertices as the size
362 const auto thisNumRows = std::max(this->numRows_, maxRowIdx + 1);
363 const auto thisNumCols = std::max(this->numCols_, maxColIdx + 1);
364
365 this->preparePushbackRowGrouping(thisNumRows, i);
366
367 this->groupAndTrackColumnIndicesByRow(i, j);
368
369 if constexpr (TrackCompressedIdx) {
370 if (expandExistingIdxMap) {
371 this->remapCompressedIndex(std::move(compressedIdx), numOrigNNZ);
372 }
373 }
374
375 this->numRows_ = thisNumRows;
376 this->numCols_ = thisNumCols;
377}
378
379template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
380void
382CSR::compress(const Offset maxNumVertices)
383{
384 if (this->numRows() > maxNumVertices) {
385 throw std::invalid_argument {
386 "Number of vertices in input graph (" +
387 std::to_string(this->numRows()) + ") "
388 "exceeds maximum graph size implied by explicit size of "
389 "adjacency matrix (" + std::to_string(maxNumVertices) + ')'
390 };
391 }
392
393 this->sortColumnIndicesPerRow();
394
395 // Must be called *after* sortColumnIndicesPerRow().
396 this->condenseDuplicates();
397
398 const auto nRows = this->startPointers().size() - 1;
399 if (nRows < maxNumVertices) {
400 this->ia_.insert(this->ia_.end(),
401 maxNumVertices - nRows,
402 this->startPointers().back());
403 }
404}
405
406template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
407void
410{
411 // Transposition is, in this context, effectively a linear time (O(nnz))
412 // bucket insertion procedure. In other words transposing the structure
413 // twice creates a structure with column indices in (ascendingly) sorted
414 // order.
415
416 this->transpose();
417 this->transpose();
418}
419
420template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
421void
424{
425 // Note: Must be called *after* sortColumnIndicesPerRow().
426
427 const auto colIdx = this->ja_;
428 auto end = colIdx.begin();
429
430 this->ja_.clear();
431
432 [[maybe_unused]] auto compressedIdx = this->compressedIdx_;
433 if constexpr (TrackCompressedIdx) {
434 this->compressedIdx_.clear();
435 }
436
437 const auto numRows = this->ia_.size() - 1;
438 for (auto row = 0*numRows; row < numRows; ++row) {
439 auto begin = end;
440
441 std::advance(end, this->ia_[row + 1] - this->ia_[row + 0]);
442
443 const auto q = this->ja_.size();
444
445 this->condenseAndTrackUniqueColumnsForSingleRow(begin, end);
446
447 this->ia_[row + 0] = q;
448 }
449
450 if constexpr (TrackCompressedIdx) {
451 this->remapCompressedIndex(std::move(compressedIdx));
452 }
453
454 // Record final table sizes.
455 this->ia_.back() = this->ja_.size();
456}
457
458template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
459void
461CSR::preparePushbackRowGrouping(const int numRows,
462 const Neighbours& rowIdx)
463{
464 assert (numRows >= 0);
465
466 this->ia_.assign(numRows + 1, 0);
467
468 // Count number of neighbouring vertices for each row. Accumulate in
469 // "next" bin since we're positioning the end pointers.
470 for (const auto& row : rowIdx) {
471 this->ia_[row + 1] += 1;
472 }
473
474 // Position "end" pointers.
475 //
476 // After this loop, ia_[i + 1] points to the *start* of the range of the
477 // column indices/neighbouring vertices of vertex 'i'. This, in turn,
478 // enables using the statement ja_[ia_[i+1]++] = v in groupAndTrack()
479 // to insert vertex 'v' as a neighbour, at the end of the range of known
480 // neighbours, *and* advance the end pointer of row/vertex 'i'. We use
481 // ia_[0] as an accumulator for the total number of neighbouring
482 // vertices in the graph.
483 //
484 // Note index range: 1..numRows inclusive.
485 for (typename Start::size_type i = 1, n = numRows; i <= n; ++i) {
486 this->ia_[0] += this->ia_[i];
487 this->ia_[i] = this->ia_[0] - this->ia_[i];
488 }
489
490 assert (this->ia_[0] == rowIdx.size());
491}
492
493template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
494void
496CSR::groupAndTrackColumnIndicesByRow(const Neighbours& rowIdx,
497 const Neighbours& colIdx)
498{
499 assert (this->ia_[0] == rowIdx.size());
500
501 const auto nnz = rowIdx.size();
502
503 this->ja_.resize(nnz);
504
505 if constexpr (TrackCompressedIdx) {
506 this->compressedIdx_.clear();
507 this->compressedIdx_.reserve(nnz);
508 }
509
510 // Group/insert column indices according to their associate vertex/row
511 // index.
512 //
513 // At the start of the loop the end pointers ia_[i+1], formed in
514 // preparePushback(), are positioned at the *start* of the column index
515 // range associated to vertex 'i'. After this loop all vertices
516 // neighbouring vertex 'i' will be placed consecutively, in order of
517 // appearance, into ja_. Furthermore, the row pointers ia_ will have
518 // their final position.
519 //
520 // The statement ja_[ia_[i+1]++] = v, split into two statements using
521 // the helper object 'k', inserts 'v' as a neighbouring vertex of vertex
522 // 'i' *and* advances the end pointer ia_[i+1] of that vertex. We use
523 // and maintain the invariant that ia_[i+1] at all times records the
524 // insertion point of the next neighbouring vertex of vertex 'i'. When
525 // the list of neighbouring vertices for vertex 'i' has been exhausted,
526 // ia_[i+1] will hold the start position for in ja_ for vertex i+1.
527 for (auto nz = 0*nnz; nz < nnz; ++nz) {
528 const auto k = this->ia_[rowIdx[nz] + 1] ++;
529
530 this->ja_[k] = colIdx[nz];
531
532 if constexpr (TrackCompressedIdx) {
533 this->compressedIdx_.push_back(k);
534 }
535 }
536
537 this->ia_[0] = 0;
538}
539
540template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
541void
544{
545 [[maybe_unused]] auto compressedIdx = this->compressedIdx_;
546
547 {
548 const auto rowIdx = this->coordinateFormatRowIndices();
549 const auto colIdx = this->ja_;
550
551 this->preparePushbackRowGrouping(this->numCols_, colIdx);
552
553 // Note parameter order. Transposition switches role of rows and
554 // columns.
555 this->groupAndTrackColumnIndicesByRow(colIdx, rowIdx);
556 }
557
558 if constexpr (TrackCompressedIdx) {
559 this->remapCompressedIndex(std::move(compressedIdx));
560 }
561
562 std::swap(this->numRows_, this->numCols_);
563}
564
565template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
566void
568CSR::condenseAndTrackUniqueColumnsForSingleRow(typename Neighbours::const_iterator begin,
569 typename Neighbours::const_iterator end)
570{
571 // We assume that we're only called *after* sortColumnIndicesPerRow()
572 // whence duplicate elements appear consecutively in [begin, end).
573 //
574 // Note: This is essentially the same as std::unique(begin, end) save
575 // for the return value and the fact that we additionally record the
576 // 'compressedIdx_' mapping. That mapping enables subsequent, decoupled
577 // accumulation of the 'sa_' contributions.
578
579 while (begin != end) {
580 // Note: Order of ja_ and compressedIdx_ matters here.
581
582 if constexpr (TrackCompressedIdx) {
583 this->compressedIdx_.push_back(this->ja_.size());
584 }
585
586 this->ja_.push_back(*begin);
587
588 auto next_unique =
589 std::find_if(begin, end, [last = this->ja_.back()]
590 (const auto j) { return j != last; });
591
592 if constexpr (TrackCompressedIdx) {
593 // Number of duplicate elements in [begin, next_unique).
594 const auto ndup = std::distance(begin, next_unique);
595
596 if (ndup > 1) {
597 // Insert ndup - 1 copies of .back() to represent the
598 // duplicate pairs in [begin, next_unique). We subtract one
599 // to account for .push_back() above representing *begin.
600 this->compressedIdx_.insert(this->compressedIdx_.end(),
601 ndup - 1,
602 this->compressedIdx_.back());
603 }
604 }
605
606 begin = next_unique;
607 }
608}
609
610template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
611void
613remapCompressedIndex([[maybe_unused]] Start&& compressedIdx,
614 [[maybe_unused]] std::optional<typename Start::size_type> numOrig)
615{
616 if constexpr (TrackCompressedIdx) {
617 std::ranges::transform(compressedIdx, compressedIdx.begin(),
618 [this](const auto& i)
619 { return this->compressedIdx_[i]; });
620
621 if (numOrig.has_value() && (*numOrig < this->compressedIdx_.size())) {
622 // Client called add() after compress(). Remap existing portion
623 // of compressedIdx (above), and append new entries (here).
624 compressedIdx
625 .insert(compressedIdx.end(),
626 this->compressedIdx_.begin() + *numOrig,
627 this->compressedIdx_.end());
628 }
629
630 this->compressedIdx_.swap(compressedIdx);
631 }
632}
633
634// =====================================================================
635
636// ---------------------------------------------------------------------
637// Class Opm::utility::CSRGraphFromCoordinates
638// ---------------------------------------------------------------------
639
640template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
646
647template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
648void
650addConnection(const VertexID v1, const VertexID v2)
651{
652 if ((v1 < 0) || (v2 < 0)) {
653 throw std::invalid_argument {
654 "Vertex IDs must be non-negative. Got (v1,v2) = ("
655 + std::to_string(v1) + ", " + std::to_string(v2)
656 + ')'
657 };
658 }
659
660 if constexpr (! PermitSelfConnections) {
661 if (v1 == v2) {
662 // Ignore self connections.
663 return;
664 }
665 }
666
667 this->uncompressed_.add(v1, v2);
668}
669
670template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
671void
673addVertexGroup(const std::vector<VertexID>& vertices)
674{
675 if (vertices.empty()) {
676 return;
677 }
678 // Initialize any new vertices in the disjoint set union structure
679 for (const auto& v : vertices) {
680 if (parent_.find(v) == parent_.end()) {
681 parent_.emplace(v, v); // Each vertex initially points to itself
682 }
683 }
684
685 // Union all vertices in the group
686 if (vertices.size() > 1) {
687 for (size_t i = 1; i < vertices.size(); ++i) {
688 unionSets(vertices[0], vertices[i]);
689 }
690 }
691}
692
693template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
694VertexID
696find(VertexID v)
697{
698 // If vertex is not in the structure, add it as its own parent
699 auto it = parent_.find(v);
700 if (it == parent_.end()) {
701 parent_.emplace(v, v);
702 return v;
703 }
704
705 // Path compression: update parent pointers to point directly to the root
706 if (it->second != v) {
707 it->second = find(it->second);
708 }
709 return it->second;
710}
711
712template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
713void
714Opm::utility::CSRGraphFromCoordinates<VertexID, TrackCompressedIdx, PermitSelfConnections>::
715unionSets(VertexID a, VertexID b)
716{
717 // Find representatives (roots) of both sets
718 VertexID rootA = find(a);
719 VertexID rootB = find(b);
720
721 // If already in same set, do nothing
722 if (rootA == rootB) {
723 return;
724 }
725
726 // Union by rank: we'll use the smaller vertex ID as the root
727 // (This is a simple heuristic that works well for this use case)
728 if (rootA < rootB) {
729 parent_.at(rootB) = rootA;
730 } else {
731 parent_.at(rootA) = rootB;
732 }
733}
734
735template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
739{
740 // If no vertex groups were defined, nothing to do
741 if (parent_.empty()) {
742 return this->uncompressed_.maxRow().value_or(0) + 1;
743 }
744
745 // Create the direct mapping for merges (original → target)
746 std::unordered_map<VertexID, VertexID> vertex_merges;
747 for (auto& [vertex, parent] : parent_) {
748 // Find the final root for this vertex which is the min vertex ID in the set
749 VertexID root = find(vertex);
750
751 // Only add to merges if the vertex isn't its own root
752 if (vertex != root) {
753 vertex_merges.emplace(vertex, root);
754 }
755 }
756 // Apply the merges to the uncompressed data
757 if (!vertex_merges.empty()) {
758 vertex_mapping_ = this->uncompressed_.applyVertexMerges(vertex_merges);
759 }
760
761 return this->uncompressed_.maxRow().value_or(0) + 1;
762}
763
764template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
765void
767compress(Offset maxNumVertices, bool expandExistingIdxMap)
768{
769 // Apply vertex merges if there are vertex groups but merges haven't been applied yet
770 if (!parent_.empty() && vertex_mapping_.empty()) {
772 }
773
774 if (! this->uncompressed_.isValid()) {
775 throw std::logic_error {
776 "Cannot compress invalid connection list"
777 };
778 }
779
780 this->csr_.merge(this->uncompressed_, maxNumVertices, expandExistingIdxMap);
781
782 this->uncompressed_.clear();
783}
784
785template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
786VertexID
788getFinalVertexID(VertexID v) const
789{
790 if (vertex_mapping_.empty()) {
791 return v;
792 }
793 else {
794 return vertex_mapping_.at(v);
795 }
796}
797
798template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
804
805template <typename VertexID, bool TrackCompressedIdx, bool PermitSelfConnections>
808{
809 const auto& ia = this->startPointers();
810
811 return ia.empty() ? 0 : ia.back();
812}
Form CSR adjacency matrix representation of unstructured graphs.
Definition CSRGraphFromCoordinates.hpp:56
const Start & startPointers() const
Read-only access to compressed structure's start pointers.
Definition CSRGraphFromCoordinates.hpp:121
Offset numEdges() const
Retrieve number of edges (non-zero matrix elements) in input graph.
Definition CSRGraphFromCoordinates_impl.hpp:807
void clear()
Clear all internal buffers, but preserve allocated capacity.
Definition CSRGraphFromCoordinates_impl.hpp:641
void addConnection(VertexID v1, VertexID v2)
Add flow rate connection between regions.
Definition CSRGraphFromCoordinates_impl.hpp:650
VertexID getFinalVertexID(VertexID originalVertexID) const
Get the final vertex ID after all merges and renumbering for a given original vertex ID.
Definition CSRGraphFromCoordinates_impl.hpp:788
typename Neighbours::size_type Offset
Offset into neighbour array.
Definition CSRGraphFromCoordinates.hpp:68
Offset numVertices() const
Retrieve number of rows (source entities) in input graph.
Definition CSRGraphFromCoordinates_impl.hpp:800
void addVertexGroup(const std::vector< VertexID > &vertices)
Add a group of vertices that should be merged together.
Definition CSRGraphFromCoordinates_impl.hpp:673
Offset applyVertexMerges()
Apply vertex merges to all vertex groups.
Definition CSRGraphFromCoordinates_impl.hpp:738
void compress(Offset maxNumVertices, bool expandExistingIdxMap=false)
Form CSR adjacency matrix representation of input graph from connections established in previous call...
Definition CSRGraphFromCoordinates_impl.hpp:767
std::vector< BaseVertexID > Neighbours
Representation of neighbouring regions.
Definition CSRGraphFromCoordinates.hpp:65
std::vector< Offset > Start
CSR start pointers.
Definition CSRGraphFromCoordinates.hpp:71