pointDataI.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) 2011-2022 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
34  s_(great),
35  v_(point::max)
36 {}
37 
38 
40 (
41  const point& origin,
42  const scalar distSqr,
43  const scalar s,
44  const vector& v
45 )
46 :
47  pointEdgePoint(origin, distSqr),
48  s_(s),
49  v_(v)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 inline Foam::scalar Foam::pointData::s() const
56 {
57  return s_;
58 }
59 
60 
61 inline const Foam::vector& Foam::pointData::v() const
62 {
63  return v_;
64 }
65 
66 
67 template<class TrackingData>
69 (
70  const polyPatch& patch,
71  const label patchFacei,
72  const transformer& transform,
73  TrackingData& td
74 )
75 {
76  pointEdgePoint::transform(patch, patchFacei, transform, td);
77  v_ = transform.transform(v_);
78 }
79 
80 
81 // Update this with information from connected edge
82 template<class TrackingData>
84 (
85  const polyMesh& mesh,
86  const label pointi,
87  const label edgei,
88  const pointData& edgeinfo,
89  const scalar tol,
90  TrackingData& td
91 )
92 {
93  if
94  (
96  (
97  mesh,
98  pointi,
99  edgei,
100  edgeinfo,
101  tol,
102  td
103  )
104  )
105  {
106  s_ = edgeinfo.s_;
107  v_ = edgeinfo.v_;
108  return true;
109  }
110  else
111  {
112  return false;
113  }
114 }
115 
116 // Update this with new information on same point
117 template<class TrackingData>
119 (
120  const polyMesh& mesh,
121  const label pointi,
122  const pointData& newPointInfo,
123  const scalar tol,
124  TrackingData& td
125 )
126 {
127  if
128  (
130  (
131  mesh,
132  pointi,
133  newPointInfo,
134  tol,
135  td
136  )
137  )
138  {
139  s_ = newPointInfo.s_;
140  v_ = newPointInfo.v_;
141  return true;
142  }
143  else
144  {
145  return false;
146  }
147 }
148 
149 
150 // Update this with new information on same point. No extra information.
151 template<class TrackingData>
153 (
154  const pointData& newPointInfo,
155  const scalar tol,
156  TrackingData& td
157 )
158 {
159  if (pointEdgePoint::updatePoint(newPointInfo, tol, td))
160  {
161  s_ = newPointInfo.s_;
162  v_ = newPointInfo.v_;
163  return true;
164  }
165  else
166  {
167  return false;
168  }
169 }
170 
171 
172 // Update this with information from connected point
173 template<class TrackingData>
175 (
176  const polyMesh& mesh,
177  const label edgei,
178  const label pointi,
179  const pointData& pointInfo,
180  const scalar tol,
181  TrackingData& td
182 
183 )
184 {
185  if
186  (
188  (
189  mesh,
190  edgei,
191  pointi,
192  pointInfo,
193  tol,
194  td
195  )
196  )
197  {
198  s_ = pointInfo.s_;
199  v_ = pointInfo.v_;
200  return true;
201  }
202  else
203  {
204  return false;
205  }
206 }
207 
208 
209 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
210 
212 const
213 {
214  return
216  && (s() == rhs.s())
217  && (v() == rhs.v());
218 }
219 
220 
222 const
223 {
224  return !(*this == rhs);
225 }
226 
227 
228 // ************************************************************************* //
friend dimensionSet transform(const dimensionSet &)
Return the argument; transformations do not change the dimensions.
Variant of pointEdgePoint with some transported additional data. WIP - should be templated on data li...
Definition: pointData.H:63
pointData()
Construct null.
Definition: pointDataI.H:31
bool updateEdge(const polyMesh &mesh, const label edgei, const label pointi, const pointData &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Definition: pointDataI.H:175
const vector & v() const
Definition: pointDataI.H:61
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgei, const pointData &edgeinfo, const scalar tol, TrackingData &td)
Influence of edge on point.
Definition: pointDataI.H:84
bool operator==(const pointData &) const
Definition: pointDataI.H:211
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Definition: pointDataI.H:69
bool operator!=(const pointData &) const
Definition: pointDataI.H:221
scalar s() const
Definition: pointDataI.H:55
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool operator==(const pointEdgePoint &) const
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
3D tensor transformation operations.