searchableSurfaceToFaceZone.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2021 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 
27 #include "polyMesh.H"
28 #include "faceZoneSet.H"
29 #include "searchableSurface.H"
30 #include "syncTools.H"
31 #include "Time.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(searchableSurfaceToFaceZone, 0);
41  (
42  topoSetSource,
43  searchableSurfaceToFaceZone,
44  word
45  );
46 }
47 
48 
49 Foam::topoSetSource::addToUsageTable Foam::searchableSurfaceToFaceZone::usage_
50 (
51  searchableSurfaceToFaceZone::typeName,
52  "\n Usage: searchableSurfaceToFaceZone surface\n\n"
53  " Select all faces whose cell-cell centre vector intersects the surface "
54  "\n"
55 );
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 void Foam::searchableSurfaceToFaceZone::combine
61 (
62  faceZoneSet& fzSet,
63  const bool add
64 ) const
65 {
66  // Get a list of "interior" faces; i.e., either internal or coupled
67  labelList interiorFaceFaces(mesh_.nFaces());
68  {
69  forAll(mesh_.faceNeighbour(), facei)
70  {
71  interiorFaceFaces[facei] = facei;
72  }
73  label nInteriorFaces = mesh_.nInternalFaces();
74  forAll(mesh_.boundaryMesh(), patchi)
75  {
76  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
77  if (patch.coupled())
78  {
79  forAll(patch, patchFacei)
80  {
81  const label facei = patch.start() + patchFacei;
82  interiorFaceFaces[nInteriorFaces] = facei;
83  ++ nInteriorFaces;
84  }
85  }
86  }
87  interiorFaceFaces.resize(nInteriorFaces);
88  }
89 
90  // Create owner and neighbour cell centre lists for all interior faces
91  pointField ownCc(interiorFaceFaces.size());
92  pointField nbrCc(interiorFaceFaces.size());
93  {
94  const pointField& cc = mesh_.cellCentres();
95 
96  vectorField boundaryNbrCc;
97  syncTools::swapBoundaryCellPositions(mesh_, cc, boundaryNbrCc);
98 
99  forAll(ownCc, interiorFacei)
100  {
101  const label facei = interiorFaceFaces[interiorFacei];
102 
103  ownCc[interiorFacei] = cc[mesh_.faceOwner()[facei]];
104  nbrCc[interiorFacei] =
105  mesh_.isInternalFace(facei)
106  ? cc[mesh_.faceNeighbour()[facei]]
107  : boundaryNbrCc[facei - mesh_.nInternalFaces()];
108  }
109  }
110 
111  // Do intersection tests on the vectors between the owner and neighbour
112  // cell centres, extended by the tolerance
113  List<pointIndexHit> hits;
114  pointField normals;
115  surfacePtr_().findLine
116  (
117  ownCc + tol_*(ownCc - nbrCc),
118  nbrCc - tol_*(ownCc - nbrCc),
119  hits
120  );
121  surfacePtr_().getNormal(hits, normals);
122 
123  // Create a list of labels indicating what side of the surface a cell
124  // is on; -1 is below, +1 is above, and 0 is too far from the surface
125  // for the sidedness to be calculable
126  labelList side(mesh_.nCells(), 0);
127  forAll(hits, interiorFacei)
128  {
129  if (hits[interiorFacei].hit())
130  {
131  const label facei = interiorFaceFaces[interiorFacei];
132 
133  const vector d = nbrCc[interiorFacei] - ownCc[interiorFacei];
134  const bool sign = (normals[interiorFacei] & d) < 0;
135 
136  side[mesh_.faceOwner()[facei]] = sign ? -1 : +1;
137 
138  if (mesh_.isInternalFace(facei))
139  {
140  side[mesh_.faceNeighbour()[facei]] = sign ? +1 : -1;
141  }
142  }
143  }
144  labelList boundaryNbrSide;
145  syncTools::swapBoundaryCellList(mesh_, side, boundaryNbrSide);
146 
147  // Create a face direction list indicating which direction a given face
148  // intersects the surface; -1 is backward, +1 is forward, and 0
149  // indicates that the face does not intersect the surface
150  labelList direction(mesh_.nFaces(), 0);
151  forAll(interiorFaceFaces, interiorFacei)
152  {
153  const label facei = interiorFaceFaces[interiorFacei];
154 
155  const label ownSide = side[mesh_.faceOwner()[facei]];
156  const label nbrSide =
157  mesh_.isInternalFace(facei)
158  ? side[mesh_.faceNeighbour()[facei]]
159  : boundaryNbrSide[facei - mesh_.nInternalFaces()];
160 
161  direction[facei] = ownSide*nbrSide < 0 ? ownSide - nbrSide : 0;
162  }
163 
164  // Select intersected faces
165  DynamicList<label> newAddressing;
166  DynamicList<bool> newFlipMap;
167  if (add)
168  {
169  // Start from copy
170  newAddressing = DynamicList<label>(fzSet.addressing());
171  newFlipMap = DynamicList<bool>(fzSet.flipMap());
172 
173  // Add everything with a direction that is not already in the set
174  forAll(direction, facei)
175  {
176  if (direction[facei] != 0 && !fzSet.found(facei))
177  {
178  newAddressing.append(facei);
179  newFlipMap.append(direction[facei] < 0);
180  }
181  }
182  }
183  else
184  {
185  // Start from empty
186  newAddressing = DynamicList<label>(fzSet.addressing().size());
187  newFlipMap = DynamicList<bool>(fzSet.flipMap().size());
188 
189  // Add everything from the zone set that does not have a direction
190  forAll(fzSet.addressing(), i)
191  {
192  if (direction[fzSet.addressing()[i]] == 0)
193  {
194  newAddressing.append(fzSet.addressing()[i]);
195  newFlipMap.append(fzSet.flipMap()[i]);
196  }
197  }
198  }
199  fzSet.addressing().transfer(newAddressing);
200  fzSet.flipMap().transfer(newFlipMap);
201  fzSet.updateSet();
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
208 (
209  const polyMesh& mesh,
210  const dictionary& dict
211 )
212 :
213  topoSetSource(mesh),
214  surfacePtr_
215  (
217  (
218  word(dict.lookup("surface")),
219  IOobject
220  (
221  dict.lookupOrDefault("name", mesh.objectRegistry::db().name()),
222  mesh.time().constant(),
224  mesh.objectRegistry::db(),
227  ),
228  dict
229  )
230  ),
231  tol_(dict.lookupOrDefault<scalar>("tol", rootSmall))
232 {}
233 
234 
235 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
236 
238 {}
239 
240 
241 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242 
244 (
245  const topoSetSource::setAction action,
246  topoSet& set
247 ) const
248 {
249  if (!isA<faceZoneSet>(set))
250  {
252  << "Operation only allowed on a faceZoneSet." << endl;
253  }
254  else
255  {
256  faceZoneSet& fzSet = refCast<faceZoneSet>(set);
257 
258  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
259  {
260  Info<< " Adding all faces from surface "
261  << surfacePtr_().name() << " ..." << endl;
262 
263  combine(fzSet, true);
264  }
265  else if (action == topoSetSource::DELETE)
266  {
267  Info<< " Removing all faces from surface "
268  << surfacePtr_().name() << " ..." << endl;
269 
270  combine(fzSet, false);
271  }
272  }
273 }
274 
275 
276 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 word & name() const
Return name.
Definition: IOobject.H:303
searchableSurfaceToFaceZone(const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static void swapBoundaryCellPositions(const polyMesh &mesh, const UList< point > &cellData, List< point > &neighbourCellData)
Swap to obtain neighbour cell positions for all boundary faces.
Definition: syncTools.C:31
AccessType combine(const List< T > &, AccessOp aop=accessOp< T >())
Combines sublists into one big list.
Definition: ListListOps.C:34
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
static const word & geometryDir()
Return the geometry directory name.
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
List< label > labelList
A List of labels.
Definition: labelList.H:56
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
Like faceSet but -reads data from faceZone -updates faceZone when writing.
Definition: faceZoneSet.H:49
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
Class with constructor to add usage string to table.
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
messageStream Info
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844