Go to the documentation of this file.
35 namespace decompositionConstraints
59 if (decompositionConstraint::debug)
61 Info<<
type() <<
" : setting constraints to preserve baffles"
73 if (decompositionConstraint::debug)
75 Info<<
type() <<
" : setting constraints to preserve baffles"
98 if (decompositionConstraint::debug & 2)
100 Info<<
type() <<
" : setting constraints to preserve "
102 <<
" baffles" <<
endl;
110 forAll(explicitConnections, i)
136 <<
"When adding baffle between faces "
137 <<
p[0] <<
" at " <<
mesh.faceCentres()[
p[0]]
139 <<
p[1] <<
" at " <<
mesh.faceCentres()[
p[1]]
140 <<
" : face " <<
p[0] <<
" already is connected to face "
141 << p0Slave <<
" at " <<
mesh.faceCentres()[p0Slave]
142 <<
" and face " <<
p[1] <<
" already is connected to face "
143 << p1Slave <<
" at " <<
mesh.faceCentres()[p1Slave]
153 if (otherFaceI != -1 && faceI < otherFaceI)
164 if (otherFaceI != -1 && faceI < otherFaceI)
166 explicitConnections[
n++] =
labelPair(faceI, otherFaceI);
173 forAll(explicitConnections, i)
175 blockedFace[explicitConnections[i].first()] =
false;
176 blockedFace[explicitConnections[i].second()] =
false;
201 const labelPair& baffle = bafflePairs[i];
205 const label procI = decomposition[
mesh.faceOwner()[f0]];
207 if (
mesh.isInternalFace(f0))
210 if (decomposition[nei0] != procI)
212 decomposition[nei0] = procI;
218 if (decomposition[own1] != procI)
220 decomposition[own1] = procI;
223 if (
mesh.isInternalFace(
f1))
226 if (decomposition[nei1] != procI)
228 decomposition[nei1] = procI;
233 if (decompositionConstraint::debug & 2)
236 Info<<
type() <<
" : changed decomposition on " << nChanged
A class for handling words, derived from string.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
addToRunTimeSelectionTable(decompositionConstraint, preserveBafflesConstraint, dictionary)
#define forAll(list, i)
Loop across all elements in list.
defineTypeName(preserveBafflesConstraint)
const Type & first() const
Return first.
virtual void add(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Add my constraints to list of constraints.
Ostream & endl(Ostream &os)
Add newline and flush stream.
preserveBafflesConstraint()
Construct from components.
Mesh consisting of general polyhedral cells.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A topoSetSource to select faces based on usage in another faceSet.
virtual void apply(const polyMesh &mesh, const boolList &blockedFace, const PtrList< labelList > &specifiedProcessorFaces, const labelList &specifiedProcessor, const List< labelPair > &explicitConnections, labelList &decomposition) const
Apply any additional post-decomposition constraints.
const Type & second() const
Return second.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
void setSize(const label)
Reset size of List.
An ordered pair of two objects of type <T> with first() and second() elements.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void size(const label)
Override size to be inconsistent with allocated storage.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Pair< label > labelPair
Label pair.
static int compare(const Pair< Type > &a, const Pair< Type > &b)
Compare Pairs.
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.