searchableSurfaceSelection.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-2025 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 "searchableSurface.H"
30 #include "fvMesh.H"
31 #include "Time.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace faceSelections
38 {
39  defineTypeNameAndDebug(searchableSurfaceSelection, 0);
41  (
42  faceSelection,
43  searchableSurfaceSelection,
44  dictionary
45  );
47  (
48  faceSelection,
49  searchableSurfaceSelection,
50  dictionary,
51  searchableSurface,
52  "searchableSurface"
53  );
54 }
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const word& name,
63  const fvMesh& mesh,
64  const dictionary& dict
65 )
66 :
67  faceSelection(name, mesh, dict),
68  surfacePtr_
69  (
70  searchableSurface::New
71  (
72  word(dict.lookup("surface")),
73  IOobject
74  (
75  dict.lookupOrDefault("name", mesh.objectRegistry::db().name()),
76  mesh.time().constant(),
77  searchableSurface::geometryDir(mesh.time()),
78  mesh.objectRegistry::db(),
79  IOobject::MUST_READ,
80  IOobject::NO_WRITE
81  ),
82  dict
83  )
84  )
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 (
98  const label zoneID,
99  labelList& faceToZoneID,
100  boolList& faceToFlip
101 ) const
102 {
103  // Get cell-cell centre vectors
104 
105  pointField start(mesh_.nFaces());
106  pointField end(mesh_.nFaces());
107 
108  // Internal faces
109  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
110  {
111  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
112  end[facei] = mesh_.cellCentres()[mesh_.faceNeighbour()[facei]];
113  }
114 
115  // Boundary faces
116  vectorField neighbourCellCentres;
118  (
119  mesh_,
120  mesh_.cellCentres(),
121  neighbourCellCentres
122  );
123 
124  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
125 
126  forAll(pbm, patchi)
127  {
128  const polyPatch& pp = pbm[patchi];
129 
130  if (pp.coupled())
131  {
132  forAll(pp, i)
133  {
134  label facei = pp.start()+i;
135  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
136  end[facei] = neighbourCellCentres[facei-mesh_.nInternalFaces()];
137  }
138  }
139  else
140  {
141  forAll(pp, i)
142  {
143  label facei = pp.start()+i;
144  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
145  end[facei] = mesh_.faceCentres()[facei];
146  }
147  }
148  }
149 
150  List<pointIndexHit> hits;
151  surfacePtr_().findLine(start, end, hits);
152  pointField normals;
153  surfacePtr_().getNormal(hits, normals);
154 
155  //- Note: do not select boundary faces.
156 
157  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
158  {
159  if (hits[facei].hit())
160  {
161  faceToZoneID[facei] = zoneID;
162  vector d = end[facei]-start[facei];
163  faceToFlip[facei] = ((normals[facei] & d) < 0);
164  }
165  }
166  forAll(pbm, patchi)
167  {
168  const polyPatch& pp = pbm[patchi];
169 
170  if (pp.coupled())
171  {
172  forAll(pp, i)
173  {
174  label facei = pp.start()+i;
175  if (hits[facei].hit())
176  {
177  faceToZoneID[facei] = zoneID;
178  vector d = end[facei]-start[facei];
179  faceToFlip[facei] = ((normals[facei] & d) < 0);
180  }
181  }
182  }
183  }
184 
185  faceSelection::select(zoneID, faceToZoneID, faceToFlip);
186 }
187 
188 
189 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
virtual void select(const label, labelList &, boolList &) const =0
searchableSurfaceSelection(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
virtual void select(const label zoneID, labelList &, boolList &) const
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
addBackwardCompatibleToRunTimeSelectionTable(fvPatchScalarField, coupledTemperatureFvPatchScalarField, patchMapper, turbulentTemperatureCoupledBaffleMixed, "compressible::turbulentTemperatureCoupledBaffleMixed")
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
defineTypeNameAndDebug(combustionModel, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict