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-2016 OpenFOAM Foundation
6  \\/ M anipulation |
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 {
38 
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
54  const dictionary& constraintsDict,
55  const word& modelType
56 )
57 :
58  decompositionConstraint(constraintsDict, typeName)
59 {
60  if (decompositionConstraint::debug)
61  {
62  Info<< type() << " : setting constraints to preserve baffles"
63  //<< returnReduce(bafflePairs.size(), sumOp<label>())
64  << endl;
65  }
66 }
67 
68 
71 :
73 {
74  if (decompositionConstraint::debug)
75  {
76  Info<< type() << " : setting constraints to preserve baffles"
77  //<< returnReduce(bafflePairs.size(), sumOp<label>())
78  << endl;
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
84 
86 (
87  const polyMesh& mesh,
88  boolList& blockedFace,
89  PtrList<labelList>& specifiedProcessorFaces,
90  labelList& specifiedProcessor,
91  List<labelPair>& explicitConnections
92 ) const
93 {
94  const labelPairList bafflePairs
95  (
97  );
98 
99  if (decompositionConstraint::debug & 2)
100  {
101  Info<< type() << " : setting constraints to preserve "
102  << returnReduce(bafflePairs.size(), sumOp<label>())
103  << " baffles" << endl;
104  }
105 
106 
107  // Merge into explicitConnections
108  {
109  // Convert into face-to-face addressing
110  labelList faceToFace(mesh.nFaces(), -1);
111  forAll(explicitConnections, i)
112  {
113  const labelPair& p = explicitConnections[i];
114  faceToFace[p[0]] = p[1];
115  faceToFace[p[1]] = p[0];
116  }
117 
118  // Merge in bafflePairs
119  forAll(bafflePairs, i)
120  {
121  const labelPair& p = bafflePairs[i];
122 
123  if (faceToFace[p[0]] == -1 && faceToFace[p[1]] == -1)
124  {
125  faceToFace[p[0]] = p[1];
126  faceToFace[p[1]] = p[0];
127  }
128  else if (labelPair::compare(p, labelPair(p[0], faceToFace[p[0]])))
129  {
130  // Connection already present
131  }
132  else
133  {
134  label p0Slave = faceToFace[p[0]];
135  label p1Slave = faceToFace[p[1]];
137  << "When adding baffle between faces "
138  << p[0] << " at " << mesh.faceCentres()[p[0]]
139  << " and "
140  << p[1] << " at " << mesh.faceCentres()[p[1]]
141  << " : face " << p[0] << " already is connected to face "
142  << p0Slave << " at " << mesh.faceCentres()[p0Slave]
143  << " and face " << p[1] << " already is connected to face "
144  << p1Slave << " at " << mesh.faceCentres()[p1Slave]
145  << endl;
146  }
147  }
148 
149  // Convert back into labelPairList
150  label n = 0;
151  forAll(faceToFace, faceI)
152  {
153  label otherFaceI = faceToFace[faceI];
154  if (otherFaceI != -1 && faceI < otherFaceI)
155  {
156  // I am master of slave
157  n++;
158  }
159  }
160  explicitConnections.setSize(n);
161  n = 0;
162  forAll(faceToFace, faceI)
163  {
164  label otherFaceI = faceToFace[faceI];
165  if (otherFaceI != -1 && faceI < otherFaceI)
166  {
167  explicitConnections[n++] = labelPair(faceI, otherFaceI);
168  }
169  }
170  }
171 
172  // Make sure blockedFace is uptodate
173  blockedFace.setSize(mesh.nFaces(), true);
174  forAll(explicitConnections, i)
175  {
176  blockedFace[explicitConnections[i].first()] = false;
177  blockedFace[explicitConnections[i].second()] = false;
178  }
179  syncTools::syncFaceList(mesh, blockedFace, andEqOp<bool>());
180 }
181 
182 
184 (
185  const polyMesh& mesh,
186  const boolList& blockedFace,
187  const PtrList<labelList>& specifiedProcessorFaces,
188  const labelList& specifiedProcessor,
189  const List<labelPair>& explicitConnections,
190  labelList& decomposition
191 ) const
192 {
193  const labelPairList bafflePairs
194  (
196  );
197 
198  label nChanged = 0;
199 
200  forAll(bafflePairs, i)
201  {
202  const labelPair& baffle = bafflePairs[i];
203  label f0 = baffle.first();
204  label f1 = baffle.second();
205 
206  const label procI = decomposition[mesh.faceOwner()[f0]];
207 
208  if (mesh.isInternalFace(f0))
209  {
210  label nei0 = mesh.faceNeighbour()[f0];
211  if (decomposition[nei0] != procI)
212  {
213  decomposition[nei0] = procI;
214  nChanged++;
215  }
216  }
217 
218  label own1 = mesh.faceOwner()[f1];
219  if (decomposition[own1] != procI)
220  {
221  decomposition[own1] = procI;
222  nChanged++;
223  }
224  if (mesh.isInternalFace(f1))
225  {
226  label nei1 = mesh.faceNeighbour()[f1];
227  if (decomposition[nei1] != procI)
228  {
229  decomposition[nei1] = procI;
230  }
231  }
232  }
233 
234  if (decompositionConstraint::debug & 2)
235  {
236  reduce(nChanged, sumOp<label>());
237  Info<< type() << " : changed decomposition on " << nChanged
238  << " cells" << endl;
239  }
240 }
241 
242 
243 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static int compare(const Pair< Type > &a, const Pair< Type > &b)
Compare Pairs.
Definition: Pair.H:141
const Type & second() const
Return second.
Definition: Pair.H:99
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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 & first() const
Return first.
Definition: Pair.H:87
A topoSetSource to select faces based on usage in another faceSet.
Definition: faceToFace.H:48
const vectorField & faceCentres() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
T & first()
Return the first element of the list.
Definition: UListI.H:114
scalar f1
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
addToRunTimeSelectionTable(decompositionConstraint, preserveBafflesConstraint, dictionary)
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
A class for handling words, derived from string.
Definition: word.H:59
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
dictionary coeffDict_
Model coefficients dictionary.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
label nFaces() const
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
void setSize(const label)
Reset size of List.
Definition: List.C:295
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
messageStream Info
label n
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
virtual void add(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Add my constraints to list of constraints.
Namespace for OpenFOAM.
defineTypeName(preserveBafflesConstraint)