wallBoundedStreamLineParticle.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "vectorFieldIOField.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
35 (
36  const trackingData& td,
37  const point& position,
38  const label celli,
39  const label facei
40 )
41 {
42  if (celli < 0)
43  {
45  << "Invalid cell (-1)" << abort(FatalError);
46  }
47 
48  bool foundU = false;
49  vector U(Zero);
50 
51  // If current position is different
52  if
53  (
56  )
57  {
58  // Store new location
60 
61  // Scalar fields
62  sampledScalars_.resize(td.vsInterp_.size());
63  forAll(td.vsInterp_, i)
64  {
65  sampledScalars_[i].append
66  (
67  td.vsInterp_[i].interpolate(position, celli, facei)
68  );
69  }
70 
71  // Vector fields
72  sampledVectors_.resize(td.vvInterp_.size());
73  forAll(td.vvInterp_, i)
74  {
75  sampledVectors_[i].append
76  (
77  td.vvInterp_[i].interpolate(position, celli, facei)
78  );
79 
80  if (td.vvInterp_.get(i) == &(td.UInterp_))
81  {
82  foundU = true;
83  U = sampledVectors_[i].last();
84  }
85  }
86  }
87 
88  if (!foundU)
89  {
90  U = td.UInterp_.interpolate(position, celli, facei);
91  }
92 
93  return U;
94 }
95 
96 
98 (
99  trackingData& td
100 )
101 {
102  vector U = interpolateFields(td, localPosition_, cell(), face());
103 
104  const scalar magU = mag(U);
105 
106  if (magU < SMALL)
107  {
108  // Stagnant particle. Might as well stop
109  lifeTime_ = 0;
110  return vector::zero;
111  }
112 
113  if (!trackForward_)
114  {
115  U = -U;
116  }
117 
118  return U/magU;
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const polyMesh& mesh,
127  const point& position,
128  const label celli,
129  const label tetFacei,
130  const label tetPti,
131  const label meshEdgeStart,
132  const label diagEdge,
133  const bool trackForward,
134  const label lifeTime
135 )
136 :
138  (
139  mesh,
140  position,
141  celli,
142  tetFacei,
143  tetPti,
144  meshEdgeStart,
145  diagEdge
146  ),
147  trackForward_(trackForward),
148  lifeTime_(lifeTime)
149 {}
150 
151 
153 (
154  const polyMesh& mesh,
155  Istream& is,
156  bool readFields,
157  bool newFormat
158 )
159 :
160  wallBoundedParticle(mesh, is, readFields, newFormat)
161 {
162  if (readFields)
163  {
164  List<scalarList> sampledScalars;
165  List<vectorList> sampledVectors;
166 
167  is >> trackForward_ >> lifeTime_
168  >> sampledPositions_ >> sampledScalars >> sampledVectors;
169 
170  sampledScalars_.resize(sampledScalars.size());
171  forAll(sampledScalars, i)
172  {
173  sampledScalars_[i].transfer(sampledScalars[i]);
174  }
175  sampledVectors_.resize(sampledVectors.size());
176  forAll(sampledVectors, i)
177  {
178  sampledVectors_[i].transfer(sampledVectors[i]);
179  }
180  }
181 
182  is.check(FUNCTION_NAME);
183 }
184 
185 
187 (
188  const wallBoundedStreamLineParticle& p
189 )
190 :
191  wallBoundedParticle(p),
192  trackForward_(p.trackForward_),
193  lifeTime_(p.lifeTime_),
194  sampledPositions_(p.sampledPositions_),
195  sampledScalars_(p.sampledScalars_),
196  sampledVectors_(p.sampledVectors_)
197 {}
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
203 (
205 )
206 {
207  if (!c.size())
208  {
209  return;
210  }
211 
213 
214  IOField<label> lifeTime
215  (
216  c.fieldIOobject("lifeTime", IOobject::MUST_READ)
217  );
218  c.checkFieldIOobject(c, lifeTime);
219 
220  vectorFieldIOField sampledPositions
221  (
222  c.fieldIOobject("sampledPositions", IOobject::MUST_READ)
223  );
224  c.checkFieldIOobject(c, sampledPositions);
225 
226  label i = 0;
227  for (wallBoundedStreamLineParticle& p : c)
228  {
229  p.lifeTime_ = lifeTime[i];
230  p.sampledPositions_.transfer(sampledPositions[i]);
231  ++i;
232  }
233 }
234 
235 
237 (
238  const Cloud<wallBoundedStreamLineParticle>& c
239 )
240 {
242 
243  const label np = c.size();
244 
245  IOField<label> lifeTime
246  (
247  c.fieldIOobject("lifeTime", IOobject::NO_READ),
248  np
249  );
250  vectorFieldIOField sampledPositions
251  (
252  c.fieldIOobject("sampledPositions", IOobject::NO_READ),
253  np
254  );
255 
256  label i = 0;
257  for (const wallBoundedStreamLineParticle& p : c)
258  {
259  lifeTime[i] = p.lifeTime_;
260  sampledPositions[i] = p.sampledPositions_;
261  ++i;
262  }
263 
264  lifeTime.write();
265  sampledPositions.write();
266 }
267 
268 
269 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
270 
271 Foam::Ostream& Foam::operator<<
272 (
273  Ostream& os,
274  const wallBoundedStreamLineParticle& p
275 )
276 {
277  os << static_cast<const wallBoundedParticle&>(p)
278  << token::SPACE << p.trackForward_
279  << token::SPACE << p.lifeTime_
280  << token::SPACE << p.sampledPositions_
281  << token::SPACE << p.sampledScalars_
282  << token::SPACE << p.sampledVectors_;
283 
285  return os;
286 }
287 
288 
289 // ************************************************************************* //
wallBoundedStreamLineParticle(const polyMesh &c, const point &position, const label celli, const label tetFacei, const label tetPti, const label meshEdgeStart, const label diagEdge, const bool trackForward, const label lifeTime)
Construct from components.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static void writeFields(const Cloud< wallBoundedStreamLineParticle > &)
Write.
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
const PtrList< interpolation< vector > > & vvInterp_
IOField< vectorField > vectorFieldIOField
IO for a Field of vectorField.
static void readFields(CloudType &)
Read.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label lifeTime_
Lifetime of particle. Particle dies when reaches 0.
dynamicFvMesh & mesh
List< DynamicList< scalar > > sampledScalars_
Sampled scalars.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
static void writeFields(const CloudType &)
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const PtrList< interpolation< scalar > > & vsInterp_
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
DynamicList< point > sampledPositions_
Sampled positions.
static void readFields(Cloud< wallBoundedStreamLineParticle > &)
Read.
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:867
U
Definition: pEqn.H:72
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Class used to pass tracking data to the trackToEdge function.
volScalarField & p
List< DynamicList< vector > > sampledVectors_
Sampled vectors.
vector interpolateFields(const trackingData &td, const point &position, const label celli, const label facei)
vector position() const
Return current particle position.
Definition: particleI.H:283
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127