preserveFaceZonesConstraint.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 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace decompositionConstraints
35 {
36  defineTypeName(preserveFaceZonesConstraint);
37 
39  (
40  decompositionConstraint,
41  preserveFaceZonesConstraint,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const dictionary& constraintsDict,
54  const word& modelType
55 )
56 :
57  decompositionConstraint(constraintsDict, typeName),
58  zones_(coeffDict_.lookup("zones"))
59 {
60  if (decompositionConstraint::debug)
61  {
62  Info<< type() << " : adding constraints to keep owner and neighbour"
63  << " of faces in zones " << zones_
64  << " on same processor" << endl;
65  }
66 }
67 
68 
71 (
72  const wordReList& zones
73 )
74 :
76  zones_(zones)
77 {
78  if (decompositionConstraint::debug)
79  {
80  Info<< type() << " : adding constraints to keep owner and neighbour"
81  << " of faces in zones " << zones_
82  << " on same processor" << endl;
83  }
84 }
85 
86 
87 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
88 
90 (
91  const polyMesh& mesh,
92  boolList& blockedFace,
93  PtrList<labelList>& specifiedProcessorFaces,
94  labelList& specifiedProcessor,
95  List<labelPair>& explicitConnections
96 ) const
97 {
98  blockedFace.setSize(mesh.nFaces(), true);
99 
100  const faceZoneMesh& fZones = mesh.faceZones();
101 
102  const labelList zoneIDs = findStrings(zones_, fZones.names());
103 
104  label nUnblocked = 0;
105 
106  forAll(zoneIDs, i)
107  {
108  const faceZone& fz = fZones[zoneIDs[i]];
109 
110  forAll(fz, i)
111  {
112  if (blockedFace[fz[i]])
113  {
114  blockedFace[fz[i]] = false;
115  nUnblocked++;
116  }
117  }
118  }
119 
120  if (decompositionConstraint::debug & 2)
121  {
122  reduce(nUnblocked, sumOp<label>());
123  Info<< type() << " : unblocked " << nUnblocked << " faces" << endl;
124  }
125 
126  syncTools::syncFaceList(mesh, blockedFace, andEqOp<bool>());
127 }
128 
129 
131 (
132  const polyMesh& mesh,
133  const boolList& blockedFace,
134  const PtrList<labelList>& specifiedProcessorFaces,
135  const labelList& specifiedProcessor,
136  const List<labelPair>& explicitConnections,
137  labelList& decomposition
138 ) const
139 {
140  // If the decomposition has not enforced the constraint do it over
141  // here.
142 
143 
144  // Synchronise decomposition on boundary
145  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
146 
147  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
148 
149  labelList destProc(mesh.nFaces()-mesh.nInternalFaces(), labelMax);
150 
151  forAll(pbm, patchi)
152  {
153  const polyPatch& pp = pbm[patchi];
154 
155  const labelUList& faceCells = pp.faceCells();
156 
157  forAll(faceCells, i)
158  {
159  label bFaceI = pp.start()+i-mesh.nInternalFaces();
160  destProc[bFaceI] = decomposition[faceCells[i]];
161  }
162  }
163 
165 
166 
167  // Override if differing
168  // ~~~~~~~~~~~~~~~~~~~~~
169 
170  const faceZoneMesh& fZones = mesh.faceZones();
171 
172  const labelList zoneIDs = findStrings(zones_, fZones.names());
173 
174  label nChanged = 0;
175 
176  forAll(zoneIDs, i)
177  {
178  const faceZone& fz = fZones[zoneIDs[i]];
179 
180  forAll(fz, i)
181  {
182  label faceI = fz[i];
183 
184  label own = mesh.faceOwner()[faceI];
185 
186  if (mesh.isInternalFace(faceI))
187  {
188  label nei = mesh.faceNeighbour()[faceI];
189  if (decomposition[own] != decomposition[nei])
190  {
191  decomposition[nei] = decomposition[own];
192  nChanged++;
193  }
194  }
195  else
196  {
197  label bFaceI = faceI-mesh.nInternalFaces();
198  if (decomposition[own] != destProc[bFaceI])
199  {
200  decomposition[own] = destProc[bFaceI];
201  nChanged++;
202  }
203  }
204  }
205  }
206 
207  if (decompositionConstraint::debug & 2)
208  {
209  reduce(nChanged, sumOp<label>());
210  Info<< type() << " : changed decomposition on " << nChanged
211  << " cells" << endl;
212  }
213 }
214 
215 
216 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
#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
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:466
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1055
label nFaces() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Macros for easy insertion into run-time selection tables.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
addToRunTimeSelectionTable(decompositionConstraint, preserveBafflesConstraint, dictionary)
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
A class for handling words, derived from string.
Definition: word.H:59
virtual void add(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Add my constraints to list of constraints.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
static const label labelMax
Definition: label.H:62
preserveFaceZonesConstraint(const dictionary &constraintsDict, const word &type)
Construct with generic dictionary with optional entry for type.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
Foam::polyBoundaryMesh.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
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:281
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
messageStream Info
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Namespace for OpenFOAM.
defineTypeName(preserveBafflesConstraint)