singleProcessorFaceSetsConstraint.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 "faceSet.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace decompositionConstraints
36 {
37  defineTypeName(singleProcessorFaceSetsConstraint);
38 
40  (
41  decompositionConstraint,
42  singleProcessorFaceSetsConstraint,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
54  const dictionary& constraintsDict,
55  const word& modelType
56 )
57 :
58  decompositionConstraint(constraintsDict, typeName),
59  setNameAndProcs_(coeffDict_.lookup("singleProcessorFaceSets"))
60 {
61  if (decompositionConstraint::debug)
62  {
63  Info<< type()
64  << " : adding constraints to keep" << endl;
65 
66  forAll(setNameAndProcs_, setI)
67  {
68  Info<< " all cells connected to faceSet "
69  << setNameAndProcs_[setI].first()
70  << " on processor " << setNameAndProcs_[setI].second() << endl;
71  }
72  }
73 }
74 
75 
78 (
79  const List<Tuple2<word, label>>& setNameAndProcs
80 )
81 :
83  setNameAndProcs_(setNameAndProcs)
84 {
85  if (decompositionConstraint::debug)
86  {
87  Info<< type()
88  << " : adding constraints to keep" << endl;
89 
90  forAll(setNameAndProcs_, setI)
91  {
92  Info<< " all cells connected to faceSet "
93  << setNameAndProcs_[setI].first()
94  << " on processor " << setNameAndProcs_[setI].second() << endl;
95  }
96  }
97 }
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
103 (
104  const polyMesh& mesh,
105  boolList& blockedFace,
106  PtrList<labelList>& specifiedProcessorFaces,
107  labelList& specifiedProcessor,
108  List<labelPair>& explicitConnections
109 ) const
110 {
111  blockedFace.setSize(mesh.nFaces(), true);
112 
113  // Mark faces already in set
114  labelList faceToSet(mesh.nFaces(), -1);
115  forAll(specifiedProcessorFaces, setI)
116  {
117  const labelList& faceLabels = specifiedProcessorFaces[setI];
118  forAll(faceLabels, i)
119  {
120  faceToSet[faceLabels[i]] = setI;
121  }
122  }
123 
124 
125  forAll(setNameAndProcs_, setI)
126  {
127  //Info<< "Keeping all cells connected to faceSet "
128  // << setNameAndProcs_[setI].first()
129  // << " on processor " << setNameAndProcs_[setI].second() << endl;
130 
131  const label destProcI = setNameAndProcs_[setI].second();
132 
133  // Read faceSet
134  const faceSet fz(mesh, setNameAndProcs_[setI].first());
135 
136  // Check that it does not overlap with existing specifiedProcessorFaces
137  labelList nMatch(specifiedProcessorFaces.size(), 0);
138  forAllConstIter(faceSet, fz, iter)
139  {
140  label setI = faceToSet[iter.key()];
141  if (setI != -1)
142  {
143  nMatch[setI]++;
144  }
145  }
146 
147 
148  // Only store if all faces are not yet in specifiedProcessorFaces
149  // (on all processors)
150  bool store = true;
151 
152  forAll(nMatch, setI)
153  {
154  if (nMatch[setI] == fz.size())
155  {
156  // full match
157  store = false;
158  break;
159  }
160  else if (nMatch[setI] > 0)
161  {
162  // partial match
163  store = false;
164  break;
165  }
166  }
167 
168  reduce(store, andOp<bool>());
169 
170 
171  if (store)
172  {
173  specifiedProcessorFaces.append(new labelList(fz.sortedToc()));
174  specifiedProcessor.append(destProcI);
175  }
176  }
177 
178 
179  // Unblock all point connected faces
180  // 1. Mark all points on specifiedProcessorFaces
181  boolList procFacePoint(mesh.nPoints(), false);
182  forAll(specifiedProcessorFaces, setI)
183  {
184  const labelList& set = specifiedProcessorFaces[setI];
185  forAll(set, fI)
186  {
187  const face& f = mesh.faces()[set[fI]];
188  forAll(f, fp)
189  {
190  procFacePoint[f[fp]] = true;
191  }
192  }
193  }
194  syncTools::syncPointList(mesh, procFacePoint, orEqOp<bool>(), false);
195 
196  // 2. Unblock all faces on procFacePoint
197 
198  label nUnblocked = 0;
199 
200  forAll(procFacePoint, pointi)
201  {
202  if (procFacePoint[pointi])
203  {
204  const labelList& pFaces = mesh.pointFaces()[pointi];
205  forAll(pFaces, i)
206  {
207  if (blockedFace[pFaces[i]])
208  {
209  blockedFace[pFaces[i]] = false;
210  nUnblocked++;
211  }
212  }
213  }
214  }
215 
216  if (decompositionConstraint::debug & 2)
217  {
218  reduce(nUnblocked, sumOp<label>());
219  Info<< type() << " : unblocked " << nUnblocked << " faces" << endl;
220  }
221 
222  syncTools::syncFaceList(mesh, blockedFace, andEqOp<bool>());
223 }
224 
225 
227 (
228  const polyMesh& mesh,
229  const boolList& blockedFace,
230  const PtrList<labelList>& specifiedProcessorFaces,
231  const labelList& specifiedProcessor,
232  const List<labelPair>& explicitConnections,
233  labelList& decomposition
234 ) const
235 {
236  // For specifiedProcessorFaces rework the cellToProc to enforce
237  // all on one processor since we can't guarantee that the input
238  // to regionSplit was a single region.
239  // E.g. faceSet 'a' with the cells split into two regions
240  // by a notch formed by two walls
241  //
242  // \ /
243  // \ /
244  // ---a----+-----a-----
245  //
246  //
247  // Note that reworking the cellToProc might make the decomposition
248  // unbalanced.
249  label nChanged = 0;
250 
251  forAll(specifiedProcessorFaces, setI)
252  {
253  const labelList& set = specifiedProcessorFaces[setI];
254 
255  // Get the processor to use for the set
256  label procI = specifiedProcessor[setI];
257  if (procI == -1)
258  {
259  // If no processor specified use the one from the
260  // 0th element
261  if (set.size())
262  {
263  procI = decomposition[mesh.faceOwner()[set[0]]];
264  }
265  reduce(procI, maxOp<label>());
266  }
267 
268  // Get all points on the sets
269  boolList procFacePoint(mesh.nPoints(), false);
270  forAll(set, fI)
271  {
272  const face& f = mesh.faces()[set[fI]];
273  forAll(f, fp)
274  {
275  procFacePoint[f[fp]] = true;
276  }
277  }
278  syncTools::syncPointList(mesh, procFacePoint, orEqOp<bool>(), false);
279 
280  // 2. Unblock all faces on procFacePoint
281  forAll(procFacePoint, pointi)
282  {
283  if (procFacePoint[pointi])
284  {
285  const labelList& pFaces = mesh.pointFaces()[pointi];
286  forAll(pFaces, i)
287  {
288  label faceI = pFaces[i];
289 
290  label own = mesh.faceOwner()[faceI];
291  if (decomposition[own] != procI)
292  {
293  decomposition[own] = procI;
294  nChanged++;
295  }
296  if (mesh.isInternalFace(faceI))
297  {
298  label nei = mesh.faceNeighbour()[faceI];
299  if (decomposition[nei] != procI)
300  {
301  decomposition[nei] = procI;
302  nChanged++;
303  }
304  }
305  }
306  }
307  }
308  }
309 
310  if (decompositionConstraint::debug & 2)
311  {
312  reduce(nChanged, sumOp<label>());
313  Info<< type() << " : changed decomposition on " << nChanged
314  << " cells" << endl;
315  }
316 }
317 
318 
319 // ************************************************************************* //
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
#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 face labels.
Definition: faceSet.H:48
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
singleProcessorFaceSetsConstraint(const dictionary &constraintsDict, const word &type)
Construct with generic dictionary with optional entry for type.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1055
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
label nFaces() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
virtual void add(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Add my constraints to list of constraints.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(decompositionConstraint, preserveBafflesConstraint, dictionary)
A class for handling words, derived from string.
Definition: word.H:59
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:184
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1036
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
labelList f(nPoints)
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
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
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
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
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
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.
messageStream Info
const labelListList & pointFaces() const
label nPoints() const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Namespace for OpenFOAM.
defineTypeName(preserveBafflesConstraint)