preserveBafflesConstraint.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM 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
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "syncTools.H"
29 #include "localPointRegion.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace decompositionConstraints
36 {
39 (
43 );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const dictionary& constraintsDict,
54  const word& modelType
55 )
56 :
57  decompositionConstraint(constraintsDict, typeName)
58 {
59  if (decompositionConstraint::debug)
60  {
61  Info<< type() << " : setting constraints to preserve baffles"
62  //<< returnReduce(bafflePairs.size(), sumOp<label>())
63  << endl;
64  }
65 }
66 
67 
70 :
72 {
73  if (decompositionConstraint::debug)
74  {
75  Info<< type() << " : setting constraints to preserve baffles"
76  //<< returnReduce(bafflePairs.size(), sumOp<label>())
77  << endl;
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
83 
85 (
86  const polyMesh& mesh,
87  boolList& blockedFace,
88  PtrList<labelList>& specifiedProcessorFaces,
89  labelList& specifiedProcessor,
90  List<labelPair>& explicitConnections
91 ) const
92 {
93  const labelPairList bafflePairs
94  (
96  );
97 
98  if (decompositionConstraint::debug & 2)
99  {
100  Info<< type() << " : setting constraints to preserve "
101  << returnReduce(bafflePairs.size(), sumOp<label>())
102  << " baffles" << endl;
103  }
104 
105 
106  // Merge into explicitConnections
107  {
108  // Convert into face-to-face addressing
109  labelList faceToFace(mesh.nFaces(), -1);
110  forAll(explicitConnections, i)
111  {
112  const labelPair& p = explicitConnections[i];
113  faceToFace[p[0]] = p[1];
114  faceToFace[p[1]] = p[0];
115  }
116 
117  // Merge in bafflePairs
118  forAll(bafflePairs, i)
119  {
120  const labelPair& p = bafflePairs[i];
121 
122  if (faceToFace[p[0]] == -1 && faceToFace[p[1]] == -1)
123  {
124  faceToFace[p[0]] = p[1];
125  faceToFace[p[1]] = p[0];
126  }
127  else if (labelPair::compare(p, labelPair(p[0], faceToFace[p[0]])))
128  {
129  // Connection already present
130  }
131  else
132  {
133  label p0Slave = faceToFace[p[0]];
134  label p1Slave = faceToFace[p[1]];
135  IOWarningInFunction(coeffDict_)
136  << "When adding baffle between faces "
137  << p[0] << " at " << mesh.faceCentres()[p[0]]
138  << " and "
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]
144  << endl;
145  }
146  }
147 
148  // Convert back into labelPairList
149  label n = 0;
150  forAll(faceToFace, faceI)
151  {
152  label otherFaceI = faceToFace[faceI];
153  if (otherFaceI != -1 && faceI < otherFaceI)
154  {
155  // I am master of slave
156  n++;
157  }
158  }
159  explicitConnections.setSize(n);
160  n = 0;
161  forAll(faceToFace, faceI)
162  {
163  label otherFaceI = faceToFace[faceI];
164  if (otherFaceI != -1 && faceI < otherFaceI)
165  {
166  explicitConnections[n++] = labelPair(faceI, otherFaceI);
167  }
168  }
169  }
170 
171  // Make sure blockedFace is uptodate
172  blockedFace.setSize(mesh.nFaces(), true);
173  forAll(explicitConnections, i)
174  {
175  blockedFace[explicitConnections[i].first()] = false;
176  blockedFace[explicitConnections[i].second()] = false;
177  }
178  syncTools::syncFaceList(mesh, blockedFace, andEqOp<bool>());
179 }
180 
181 
183 (
184  const polyMesh& mesh,
185  const boolList& blockedFace,
186  const PtrList<labelList>& specifiedProcessorFaces,
187  const labelList& specifiedProcessor,
188  const List<labelPair>& explicitConnections,
189  labelList& decomposition
190 ) const
191 {
192  const labelPairList bafflePairs
193  (
195  );
196 
197  label nChanged = 0;
198 
199  forAll(bafflePairs, i)
200  {
201  const labelPair& baffle = bafflePairs[i];
202  label f0 = baffle.first();
203  label f1 = baffle.second();
204 
205  const label procI = decomposition[mesh.faceOwner()[f0]];
206 
207  if (mesh.isInternalFace(f0))
208  {
209  label nei0 = mesh.faceNeighbour()[f0];
210  if (decomposition[nei0] != procI)
211  {
212  decomposition[nei0] = procI;
213  nChanged++;
214  }
215  }
216 
217  label own1 = mesh.faceOwner()[f1];
218  if (decomposition[own1] != procI)
219  {
220  decomposition[own1] = procI;
221  nChanged++;
222  }
223  if (mesh.isInternalFace(f1))
224  {
225  label nei1 = mesh.faceNeighbour()[f1];
226  if (decomposition[nei1] != procI)
227  {
228  decomposition[nei1] = procI;
229  }
230  }
231  }
232 
233  if (decompositionConstraint::debug & 2)
234  {
235  reduce(nChanged, sumOp<label>());
236  Info<< type() << " : changed decomposition on " << nChanged
237  << " cells" << endl;
238  }
239 }
240 
241 
242 // ************************************************************************* //
Foam::decompositionConstraint
Definition: decompositionConstraint.H:54
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::decompositionConstraints::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionConstraint, preserveBafflesConstraint, dictionary)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::decompositionConstraints::defineTypeName
defineTypeName(preserveBafflesConstraint)
localPointRegion.H
Foam::Pair::first
const Type & first() const
Return first.
Definition: Pair.H:87
Foam::decompositionConstraints::preserveBafflesConstraint::add
virtual void add(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Add my constraints to list of constraints.
Definition: preserveBafflesConstraint.C:85
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::decompositionConstraints::preserveBafflesConstraint::preserveBafflesConstraint
preserveBafflesConstraint()
Construct from components.
Definition: preserveBafflesConstraint.C:69
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::faceToFace
A topoSetSource to select faces based on usage in another faceSet.
Definition: faceToFace.H:48
Foam::decompositionConstraints::preserveBafflesConstraint::apply
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.
Definition: preserveBafflesConstraint.C:183
Foam::Info
messageStream Info
f1
scalar f1
Definition: createFields.H:28
Foam::Pair::second
const Type & second() const
Return second.
Definition: Pair.H:99
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
preserveBafflesConstraint.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::sumOp
Definition: ops.H:162
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::decompositionConstraints::preserveBafflesConstraint
Definition: preserveBafflesConstraint.H:51
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:271
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::labelPair
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Foam::Pair::compare
static int compare(const Pair< Type > &a, const Pair< Type > &b)
Compare Pairs.
Definition: Pair.H:141
Foam::localPointRegion::findDuplicateFacePairs
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Definition: localPointRegion.C:620
Foam::andEqOp
Definition: ops.H:81