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-2013 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 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
searchableSurface.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
syncTools.H
constant
Constant dispersed-phase particle diameter model.
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:49
Foam::label
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
Foam::faceSelections::searchableSurfaceSelection::~searchableSurfaceSelection
virtual ~searchableSurfaceSelection()
Destructor.
Foam::syncTools::swapBoundaryCellPositions
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::faceSelections::searchableSurfaceSelection::select
virtual void select(const label zoneID, labelList &, boolList &) const
dict
dictionary dict
Definition: searchingEngine.H:14
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::boolList
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::faceSelection::select
virtual void select(const label, labelList &, boolList &) const =0
Foam::faceSelections::searchableSurfaceSelection::searchableSurfaceSelection
searchableSurfaceSelection(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
searchableSurfaceSelection.H
List
Definition: Test.C:19
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress