All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
patchEdgeFaceRegionI.H
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-2018 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 
26 #include "polyMesh.H"
27 #include "transform.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Update this with w2 if w2 nearer to pt.
32 template<class TrackingData>
33 inline bool Foam::patchEdgeFaceRegion::update
34 (
35  const patchEdgeFaceRegion& w2,
36  const scalar tol,
37  TrackingData& td
38 )
39 {
40  if (!w2.valid(td))
41  {
43  << "problem." << abort(FatalError);
44  }
45 
46  if (w2.region_ == -2 || region_ == -2)
47  {
48  // Blocked edge/face
49  return false;
50  }
51 
52  if (!valid(td))
53  {
54  // current not yet set so use any value
55  operator=(w2);
56  return true;
57  }
58  else
59  {
60  if (w2.region_ < region_)
61  {
62  operator=(w2);
63  return true;
64  }
65  else
66  {
67  return false;
68  }
69  }
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
75 // Null constructor
77 :
78  region_(-1)
79 {}
80 
81 
82 // Construct from origin, distance
84 (
85  const label region
86 )
87 :
88  region_(region)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  return region_;
97 }
98 
99 
100 template<class TrackingData>
101 inline bool Foam::patchEdgeFaceRegion::valid(TrackingData& td) const
102 {
103  return region_ != -1;
104 }
105 
106 
107 template<class TrackingData>
109 (
110  const polyMesh& mesh,
111  const indirectPrimitivePatch& patch,
112  const tensor& rotTensor,
113  const scalar tol,
114  TrackingData& td
115 )
116 {}
117 
118 
119 template<class TrackingData>
121 (
122  const polyMesh& mesh,
123  const indirectPrimitivePatch& patch,
124  const label edgeI,
125  const label facei,
126  const patchEdgeFaceRegion& faceInfo,
127  const scalar tol,
128  TrackingData& td
129 )
130 {
131  return update(faceInfo, tol, td);
132 }
133 
134 
135 template<class TrackingData>
137 (
138  const polyMesh& mesh,
139  const indirectPrimitivePatch& patch,
140  const patchEdgeFaceRegion& edgeInfo,
141  const bool sameOrientation,
142  const scalar tol,
143  TrackingData& td
144 )
145 {
146  return update(edgeInfo, tol, td);
147 }
148 
149 
150 template<class TrackingData>
152 (
153  const polyMesh& mesh,
154  const indirectPrimitivePatch& patch,
155  const label facei,
156  const label edgeI,
157  const patchEdgeFaceRegion& edgeInfo,
158  const scalar tol,
159  TrackingData& td
160 )
161 {
162  return update(edgeInfo, tol, td);
163 }
164 
165 
166 template<class TrackingData>
168 (
169  const patchEdgeFaceRegion& rhs,
170  TrackingData& td
171 ) const
172 {
173  return operator==(rhs);
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
178 
179 inline bool Foam::patchEdgeFaceRegion::operator==
180 (
181  const Foam::patchEdgeFaceRegion& rhs
182 ) const
183 {
184  return region() == rhs.region();
185 }
186 
187 
188 inline bool Foam::patchEdgeFaceRegion::operator!=
189 (
190  const Foam::patchEdgeFaceRegion& rhs
191 ) const
192 {
193  return !(*this == rhs);
194 }
195 
196 
197 // ************************************************************************* //
bool updateEdge(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label edgeI, const label facei, const patchEdgeFaceRegion &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
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
bool operator==(const patchEdgeFaceRegion &) const
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool updateFace(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label facei, const label edgeI, const patchEdgeFaceRegion &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on face.
patchEdgeFaceRegion()
Construct null.
void transform(const polyMesh &mesh, const indirectPrimitivePatch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
A list of faces which address into the list of points.
3D tensor transformation operations.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Transport of region for use in PatchEdgeFaceWave.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool equal(const patchEdgeFaceRegion &, TrackingData &) const
Same (like operator==)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74