1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | cfMesh: A library for mesh generation
5 \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6 \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of cfMesh.
11 cfMesh is free software; you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation; either version 3 of the License, or (at your
14 option) any later version.
16 cfMesh is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with cfMesh. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 inline void Foam::VRWGraph::checkIndex(const label i, const label j) const
28 if( (i < 0) || (i >= rows_.size()) )
32 "void Foam::VRWGraph<T,width>::"
33 "checkIndex(const label i, const label j) const"
34 ) << "Row index " << i
35 << " is not in range " << 0
36 << " and " << rows_.size() << abort(FatalError);
39 if( (j < 0) || (j >= rows_[i].size()) )
42 "void Foam::VRWGraph<T,width>::"
43 "checkIndex(label const, const label) const"
44 ) << "Column index " << j
45 << " is not in range " << 0
46 << " and " << rows_[i].size() << abort(FatalError);
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 inline Foam::VRWGraph::VRWGraph()
58 //- Construct given size
59 inline Foam::VRWGraph::VRWGraph
67 for(label rowI=0;rowI<size;++rowI)
69 rows_[rowI].start() = INVALIDROW;
70 rows_[rowI].size() = NONE;
74 inline Foam::VRWGraph::VRWGraph
77 const label nColumnsInRow
80 data_(nRows * nColumnsInRow),
83 for(label rowI=0;rowI<nRows;++rowI)
85 rows_[rowI].start() = rowI * nColumnsInRow;
86 rows_[rowI].size() = nColumnsInRow;
90 inline Foam::VRWGraph::VRWGraph
93 const label nColumnsInRow,
97 data_(nRows * nColumnsInRow, t),
100 for(label rowI=0;rowI<nRows;++rowI)
102 rows_[rowI].start() = rowI * nColumnsInRow;
103 rows_[rowI].size() = nColumnsInRow;
107 inline Foam::VRWGraph::VRWGraph
116 inline Foam::VRWGraph::~VRWGraph()
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 inline Foam::label Foam::VRWGraph::size() const
127 inline Foam::label Foam::VRWGraph::sizeOfRow(const label rowI) const
129 return rows_[rowI].size();
132 inline void Foam::VRWGraph::setSize(const label size)
134 if( size > rows_.size() )
136 rowElement rowInfo(INVALIDROW, NONE);
138 for(label i=rows_.size();i<size;++i)
139 rows_.append(rowInfo);
147 void Foam::VRWGraph::setSizeAndColumnWidth
149 const label newNumRows,
153 if( rows_.size() != 0 )
156 "void Foam::VRWGraph::setSizeAndColumnWidth"
157 "(const label size, const label rcWidth)"
158 ) << "This function should be used for empty graphs, only!"
161 data_.setSize(newNumRows * rcWidth);
164 rows_.setSize(newNumRows);
167 for(label i=0;i<newNumRows;++i)
169 rows_[i].start() = start;
171 data_[start] = FREESTART;
177 template<class ListType>
178 inline void Foam::VRWGraph::setSizeAndRowSize(const ListType& l)
180 //- set the size of graph rows
181 const label nRows = l.size();
182 rows_.setSize(nRows);
185 for(label rowI=0;rowI<nRows;++rowI)
187 rows_[rowI].size() = l[rowI];
189 if( rows_[rowI].size() != NONE )
191 rows_[rowI].start() = start;
195 rows_[rowI].start() = INVALIDROW;
198 start += rows_[rowI].size();
201 data_.setSize(start);
204 inline void Foam::VRWGraph::setRowSize(const label rowI, const label newSize)
207 if( (rowI < 0) || (rowI >= rows_.size()) )
210 "void Foam::VRWGraph<T,width>::"
211 "checkIndex(const label rowI, const label size) const"
212 ) << "Row index " << Foam::label(rowI)
213 << " is not in range " << Foam::label(0)
214 << " and " << rows_.size() << abort(FatalError);
217 const label start = rows_[rowI].start();
218 if( start == INVALIDROW )
222 rows_[rowI].start() = data_.size();
223 for(label i=0;i<newSize;++i)
225 rows_[rowI].size() = newSize;
228 else if( newSize > rows_[rowI].size() )
230 //- check if there is some unused space after the last element
231 bool foundUnused(true);
233 for(label i=rows_[rowI].size();i<newSize;++i)
235 const label j = start + i;
237 (j >= data_.size()) ||
238 (data_[j] != FREEENTRY) ||
239 (data_[j] == FREESTART)
249 //- row can be extended without copying
250 for(label i=rows_[rowI].size();i<newSize;++i)
251 data_[start+i] = NONE;
255 //- row is copied at the end of the data list
256 rows_[rowI].start() = data_.size();
257 for(label i=0;i<rows_[rowI].size();++i)
259 data_.append(data_[start+i]);
260 data_[start+i] = FREEENTRY;
262 for(label i=rows_[rowI].size();i<newSize;++i)
266 rows_[rowI].size() = newSize;
268 else if( newSize < rows_[rowI].size() )
270 for(label i=newSize;i<rows_[rowI].size();++i)
271 data_[start+i] = FREEENTRY;
272 rows_[rowI].size() = newSize;
274 rows_[rowI].start() = INVALIDROW;
278 inline void Foam::VRWGraph::clear()
284 template<class ListType>
285 inline void Foam::VRWGraph::appendList
292 rows_.append(rowElement(INVALIDROW, 0));
296 rowElement rowInfo(data_.size(), l.size());
297 const label size = l.size();
298 for(label elI=0;elI<size;++elI)
299 data_.append(l[elI]);
300 rows_.append(rowInfo);
303 inline void Foam::VRWGraph::append(const label rowI, const label el)
305 rowElement& re = rows_[rowI];
307 if( re.start() == INVALIDROW )
309 re.start() = data_.size();
315 const label oldStart = re.start();
316 const label oldSize = re.size();
319 if( oldStart + oldSize < data_.size() )
322 (data_[oldStart+oldSize] == FREEENTRY) ||
323 (data_[oldStart+oldSize] == FREESTART)
326 data_[oldStart + oldSize] = el;
330 re.start() = data_.size();
331 for(label i=0;i<oldSize;++i)
333 data_.append(data_[oldStart+i]);
334 data_[oldStart+i] = FREEENTRY;
346 inline void Foam::VRWGraph::appendIfNotIn(const label rowI, const label el)
348 if( !contains(rowI, el) )
352 template<class ListType>
353 inline void Foam::VRWGraph::setRow
359 this->setRowSize(rowI, l.size());
360 const label start = rows_[rowI].start();
361 const label size = l.size();
362 for(label elI=0;elI<size;++elI)
363 data_[start+elI] = l[elI];
366 inline void Foam::VRWGraph::mergeGraphs(const List<VRWGraph>& graphParts)
368 const label nGraphs = graphParts.size();
369 const label nRows = graphParts[0].size();
370 forAll(graphParts, i)
372 if( nRows != graphParts[i].size() )
375 "inline void Foam::VRWGraph::mergeGraphs(const List<VRWGraph>&)"
376 ) << "Cannot merge graphs" << abort(FatalError);
379 //- find the number of elements in each row
380 labelLongList nElmtsInRow(nRows);
381 for(label rowI=0;rowI<nRows;++rowI)
384 for(label i=0;i<nGraphs;++i)
385 sum += graphParts[i].sizeOfRow(rowI);
387 nElmtsInRow[rowI] = sum;
390 setSizeAndRowSize(nElmtsInRow);
392 //- Finally, assemble the merged graph
393 for(label rowI=0;rowI<nRows;++rowI)
395 forAll(graphParts, i)
397 const VRWGraph& gp = graphParts[i];
398 for(label j=0;j<gp.sizeOfRow(rowI);++j)
399 this->operator()(rowI, --nElmtsInRow[rowI]) = gp(rowI, j);
404 template<class GraphType>
405 inline void Foam::VRWGraph::reverseAddressing
408 const GraphType& origGraph
411 const label origSize = origGraph.size();
412 labelLongList nElmtsInRow(nRows);
414 for(label rowI=0;rowI<nRows;++rowI)
415 nElmtsInRow[rowI] = 0;
417 for(label rowI=0;rowI<origSize;++rowI)
419 const label rowSize = origGraph[rowI].size();
421 for(label i=0;i<rowSize;++i)
422 ++nElmtsInRow[origGraph[rowI][i]];
425 setSizeAndRowSize(nElmtsInRow);
428 //- finally fill in the data
429 for(label rowI=0;rowI<origSize;++rowI)
431 const label rowSize = origGraph[rowI].size();
433 for(label i=0;i<rowSize;++i)
435 const label el = origGraph[rowI][i];
436 this->operator()(el, nElmtsInRow[el]++) = rowI;
441 template<class GraphType>
442 inline void Foam::VRWGraph::reverseAddressing(const GraphType& origGraph)
444 const label size = origGraph.size();
447 for(label rowI=0;rowI<size;++rowI)
449 const label rowSize = origGraph[rowI].size();
450 for(label i=0;i<rowSize;++i)
451 maxValue = Foam::max(maxValue, origGraph[rowI][i]);
455 reverseAddressing(maxValue, origGraph);
458 inline void Foam::VRWGraph::reverseAddressing
461 const VRWGraph& origGraph
464 const label origSize = origGraph.size();
465 labelLongList nElmtsInRow(nRows);
467 for(label rowI=0;rowI<nRows;++rowI)
468 nElmtsInRow[rowI] = 0;
470 for(label rowI=0;rowI<origSize;++rowI)
472 const label rowSize = origGraph.sizeOfRow(rowI);
474 for(label i=0;i<rowSize;++i)
475 ++nElmtsInRow[origGraph(rowI, i)];
478 setSizeAndRowSize(nElmtsInRow);
481 //- finally fill in the data
482 for(label rowI=0;rowI<origSize;++rowI)
484 const label rowSize = origGraph.sizeOfRow(rowI);
486 for(label i=0;i<rowSize;++i)
488 const label el = origGraph(rowI, i);
489 this->operator()(el, nElmtsInRow[el]++) = rowI;
494 inline void Foam::VRWGraph::reverseAddressing(const VRWGraph& origGraph)
496 const label size = origGraph.size();
499 for(label rowI=0;rowI<size;++rowI)
501 const label rowSize = origGraph.sizeOfRow(rowI);
502 for(label i=0;i<rowSize;++i)
503 maxValue = Foam::max(maxValue, origGraph(rowI, i));
507 reverseAddressing(maxValue, origGraph);
510 inline bool Foam::VRWGraph::contains
516 const label start = rows_[rowI].start();
517 if( start == INVALIDROW )
519 const label size = rows_[rowI].size();
521 for(register label i=0;i<size;++i)
522 if( data_[start+i] == e )
528 inline Foam::label Foam::VRWGraph::containsAtPosition
534 const label start = rows_[rowI].start();
535 if( start == INVALIDROW )
538 const label size = rows_[rowI].size();
539 for(register label i=0;i<size;++i)
540 if( data_[start+i] == e )
546 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
548 inline Foam::label Foam::VRWGraph::operator()
558 return data_[rows_[i].start() + j];
562 inline Foam::label& Foam::VRWGraph::operator()
564 const label i, const label j
571 return data_[rows_[i].start() + j];
574 inline Foam::constRow Foam::VRWGraph::operator[](const label i) const
576 return constRow(*this, i);
579 inline Foam::row Foam::VRWGraph::operator[](const label i)
581 return row(*this, i);
584 inline void Foam::VRWGraph::operator=
594 // ************************************************************************* //