searchableSurfaceSelection.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) 2012-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 "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  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& name,
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
59  faceSelection(name, mesh, dict),
60  surfacePtr_
61  (
62  searchableSurface::New
63  (
64  word(dict.lookup("surface")),
65  IOobject
66  (
67  dict.lookupOrDefault("name", mesh.objectRegistry::db().name()),
68  mesh.time().constant(),
69  "triSurface",
70  mesh.objectRegistry::db(),
71  IOobject::MUST_READ,
72  IOobject::NO_WRITE
73  ),
74  dict
75  )
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 (
90  const label zoneID,
91  labelList& faceToZoneID,
92  boolList& faceToFlip
93 ) const
94 {
95  // Get cell-cell centre vectors
96 
97  pointField start(mesh_.nFaces());
98  pointField end(mesh_.nFaces());
99 
100  // Internal faces
101  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
102  {
103  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
104  end[facei] = mesh_.cellCentres()[mesh_.faceNeighbour()[facei]];
105  }
106 
107  // Boundary faces
108  vectorField neighbourCellCentres;
110  (
111  mesh_,
112  mesh_.cellCentres(),
113  neighbourCellCentres
114  );
115 
116  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
117 
118  forAll(pbm, patchi)
119  {
120  const polyPatch& pp = pbm[patchi];
121 
122  if (pp.coupled())
123  {
124  forAll(pp, i)
125  {
126  label facei = pp.start()+i;
127  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
128  end[facei] = neighbourCellCentres[facei-mesh_.nInternalFaces()];
129  }
130  }
131  else
132  {
133  forAll(pp, i)
134  {
135  label facei = pp.start()+i;
136  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
137  end[facei] = mesh_.faceCentres()[facei];
138  }
139  }
140  }
141 
142  List<pointIndexHit> hits;
143  surfacePtr_().findLine(start, end, hits);
144  pointField normals;
145  surfacePtr_().getNormal(hits, normals);
146 
147  //- Note: do not select boundary faces.
148 
149  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
150  {
151  if (hits[facei].hit())
152  {
153  faceToZoneID[facei] = zoneID;
154  vector d = end[facei]-start[facei];
155  faceToFlip[facei] = ((normals[facei] & d) < 0);
156  }
157  }
158  forAll(pbm, patchi)
159  {
160  const polyPatch& pp = pbm[patchi];
161 
162  if (pp.coupled())
163  {
164  forAll(pp, i)
165  {
166  label facei = pp.start()+i;
167  if (hits[facei].hit())
168  {
169  faceToZoneID[facei] = zoneID;
170  vector d = end[facei]-start[facei];
171  faceToFlip[facei] = ((normals[facei] & d) < 0);
172  }
173  }
174  }
175  }
176 
177  faceSelection::select(zoneID, faceToZoneID, faceToFlip);
178 }
179 
180 
181 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
dictionary dict
#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 fvMesh & mesh_
Reference to mesh.
Definition: faceSelection.H:67
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1055
label nFaces() const
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
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
List< label > labelList
A List of labels.
Definition: labelList.H:56
const vectorField & cellCentres() const
virtual void select(const label zoneID, labelList &, boolList &) const
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const vectorField & faceCentres() const
label patchi
virtual void select(const label, labelList &, boolList &) const =0
Field< vector > vectorField
Specialisation of Field<T> for vector.
searchableSurfaceSelection(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
Namespace for OpenFOAM.