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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2016 OpenFOAM Foundation
9 Copyright (C) 2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "syncTools.H"
32#include "localPointRegion.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace decompositionConstraints
39{
41
43 (
47 );
48}
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict
57)
58:
60{
61 if (decompositionConstraint::debug)
62 {
63 Info<< type()
64 << " : setting constraints to preserve baffles"
65 //<< returnReduce(bafflePairs.size(), sumOp<label>())
66 << endl;
67 }
68}
69
70
72:
74{
75 if (decompositionConstraint::debug)
76 {
77 Info<< type()
78 << " : setting constraints to preserve baffles"
79 //<< returnReduce(bafflePairs.size(), sumOp<label>())
80 << endl;
81 }
82}
83
84
85// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
86
88(
89 const polyMesh& mesh,
90 boolList& blockedFace,
91 PtrList<labelList>& specifiedProcessorFaces,
92 labelList& specifiedProcessor,
93 List<labelPair>& explicitConnections
94) const
95{
96 const labelPairList bafflePairs
97 (
99 );
100
101 if (decompositionConstraint::debug & 2)
102 {
103 Info<< type() << " : setting constraints to preserve "
104 << returnReduce(bafflePairs.size(), sumOp<label>())
105 << " baffles" << endl;
106 }
107
108
109 // Merge into explicitConnections
110 {
111 // Convert into face-to-face addressing
112 labelList faceToFace(mesh.nFaces(), -1);
113 for (const labelPair& p : explicitConnections)
114 {
115 faceToFace[p[0]] = p[1];
116 faceToFace[p[1]] = p[0];
117 }
118
119 // Merge in bafflePairs
120 for (const labelPair& p : bafflePairs)
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 const label p0Slave = faceToFace[p[0]];
134 const 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 for (const labelPair& p : explicitConnections)
174 {
175 blockedFace[p.first()] = false;
176 blockedFace[p.second()] = false;
177 }
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 for (const labelPair& baffle : bafflePairs)
200 {
201 const label f0 = baffle.first();
202 const label f1 = baffle.second();
203
204 const label procI = decomposition[mesh.faceOwner()[f0]];
205
206 if (mesh.isInternalFace(f0))
207 {
208 const label nei0 = mesh.faceNeighbour()[f0];
209 if (decomposition[nei0] != procI)
210 {
211 decomposition[nei0] = procI;
212 ++nChanged;
213 }
214 }
215
216 const label own1 = mesh.faceOwner()[f1];
217 if (decomposition[own1] != procI)
218 {
219 decomposition[own1] = procI;
220 ++nChanged;
221 }
222 if (mesh.isInternalFace(f1))
223 {
224 const label nei1 = mesh.faceNeighbour()[f1];
225 if (decomposition[nei1] != procI)
226 {
227 decomposition[nei1] = procI;
228 }
229 }
230 }
231
232 if (decompositionConstraint::debug & 2)
233 {
234 reduce(nChanged, sumOp<label>());
235 Info<< type() << " : changed decomposition on " << nChanged
236 << " cells" << endl;
237 }
238}
239
240
241// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
static int compare(const Pair< label > &a, const Pair< label > &b)
Compare Pairs.
Definition: PairI.H:31
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void apply()=0
Apply bins.
Abstract class for handling decomposition constraints.
Detects baffles and keeps owner and neighbour on same processor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A topoSetFaceSource to select faces based on usage in another faceSet.
Definition: faceToFace.H:162
Sums a given list of (at least two or more) fields and outputs the result into a new field,...
Definition: add.H:161
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
#define defineTypeName(Type)
Define the typeName.
Definition: className.H:96
volScalarField & p
dynamicFvMesh & mesh
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:57
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333