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-2023 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 template<class TrackingData>
83 (
84  const polyMesh& mesh,
85  const label pointi,
86  const label edgei,
87  const pointData& edgeinfo,
88  const scalar tol,
89  TrackingData& td
90 )
91 {
92  if
93  (
95  (
96  mesh,
97  pointi,
98  edgei,
99  edgeinfo,
100  tol,
101  td
102  )
103  )
104  {
105  s_ = edgeinfo.s_;
106  v_ = edgeinfo.v_;
107  return true;
108  }
109  else
110  {
111  return false;
112  }
113 }
114 
115 
116 template<class TrackingData>
118 (
119  const polyMesh& mesh,
120  const label pointi,
121  const pointData& newPointInfo,
122  const scalar tol,
123  TrackingData& td
124 )
125 {
126  if
127  (
129  (
130  mesh,
131  pointi,
132  newPointInfo,
133  tol,
134  td
135  )
136  )
137  {
138  s_ = newPointInfo.s_;
139  v_ = newPointInfo.v_;
140  return true;
141  }
142  else
143  {
144  return false;
145  }
146 }
147 
148 
149 template<class TrackingData>
151 (
152  const polyMesh& mesh,
153  const label edgei,
154  const label pointi,
155  const pointData& pointInfo,
156  const scalar tol,
157  TrackingData& td
158 
159 )
160 {
161  if
162  (
164  (
165  mesh,
166  edgei,
167  pointi,
168  pointInfo,
169  tol,
170  td
171  )
172  )
173  {
174  s_ = pointInfo.s_;
175  v_ = pointInfo.v_;
176  return true;
177  }
178  else
179  {
180  return false;
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
186 
188 const
189 {
190  return
192  && (s() == rhs.s())
193  && (v() == rhs.v());
194 }
195 
196 
198 const
199 {
200  return !(*this == rhs);
201 }
202 
203 
204 // ************************************************************************* //
friend dimensionSet transform(const dimensionSet &)
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:151
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:83
bool operator==(const pointData &) const
Definition: pointDataI.H:187
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:197
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:504
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
3D tensor transformation operations.