DTRMParticle.H
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) 2017-2019 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::DTRMParticle
28
29Description
30 Discrete Transfer Radiation Model (DTRM) particle
31
32SourceFiles
33 DTRMParticle.H
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef DTRMParticle_H
38#define DTRMParticle_H
39
40#include "particle.H"
41#include "IOstream.H"
42#include "autoPtr.H"
43#include "interpolationCell.H"
44#include "volFieldsFwd.H"
45#include "reflectionModel.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
51{
52
53class DTRMParticle;
54using namespace Foam::radiation;
55
56// Forward declaration of friend functions
58
59/*---------------------------------------------------------------------------*\
60 Class DTRMParticle Declaration
61\*---------------------------------------------------------------------------*/
63class DTRMParticle
64:
65 public particle
66{
67 // Private data
68
69 //- Initial position
70 point p0_;
71
72 //- Target position
73 point p1_;
74
75 //- Initial radiation intensity [W/m2]
76 scalar I0_;
77
78 //- Radiation intensity [W/m2]
79 scalar I_;
80
81 //- Area of radiation
82 scalar dA_;
83
84 //- Trasnmissive index
85 label transmissiveId_;
86
87
88public:
89
90 friend class Cloud<DTRMParticle>;
91
92 //- Class used to pass tracking data to the trackToFace function
93 class trackingData
94 :
96 {
97 // Interpolators for continuous phase fields
98
99 const interpolationCell<scalar>& aInterp_;
100 const interpolationCell<scalar>& eInterp_;
101 const interpolationCell<scalar>& EInterp_;
102 const interpolationCell<scalar>& TInterp_;
103 const interpolationCellPoint<vector>& nHatInterp_;
104
105 //- Reflected cells
106 const labelField& relfectedCells_;
107
108 //- Ptr to reflectiom model
109 UPtrList<reflectionModel> reflection_;
110
111 //- Heat source term
112 volScalarField& Q_;
113
114
115 public:
116
117 // Constructors
118
119 inline trackingData
120 (
127 const labelField&,
130 );
131
132
133 // Member functions
134
135 inline const interpolationCell<scalar>& aInterp() const;
136 inline const interpolationCell<scalar>& eInterp() const;
137 inline const interpolationCell<scalar>& EInterp() const;
138 inline const interpolationCell<scalar>& TInterp() const;
139 inline const interpolationCellPoint<vector>& nHatInterp() const;
140 inline const labelField& relfectedCells() const;
141 inline const UPtrList<reflectionModel>& reflection() const;
142
143 inline scalar& Q(label celli);
144 };
145
146
147 // Static Data Members
148
149 //- Size in bytes of the fields
150 static const std::size_t sizeofFields_;
151
152
153 //- String representation of properties
155 (
156 particle,
157 " p0"
158 + " p1"
159 + " I0"
160 + " I"
161 + " dA"
162 + " transmissiveId";
163 );
164
165
166 // Constructors
167
168 //- Construct from components, with searching for tetFace and
169 // tetPt unless disabled by doCellFacePt = false.
171 (
172 const polyMesh& mesh,
173 const vector& position,
174 const vector& targetPosition,
175 const scalar I,
176 const label cellI,
177 const scalar dA,
178 const label transmissiveId
179 );
180
181 //- Construct from components
183 (
184 const polyMesh& mesh,
186 const label celli,
187 const label tetFacei,
188 const label tetPti,
189 const vector& position,
190 const vector& targetPosition,
191 const scalar I,
192 const scalar dA,
193 const label transmissiveId
194 );
195
196 //- Construct from Istream
198 (
199 const polyMesh& mesh,
200 Istream& is,
201 bool readFields = true,
202 bool newFormat = true
203 );
204
205 //- Construct as copy
207
208
209 //- Factory class to read-construct particles used for
210 // parallel transfer
211 class iNew
212 {
213 const polyMesh& mesh_;
214
215 public:
217 iNew(const polyMesh& mesh)
218 :
219 mesh_(mesh)
220 {}
223 {
225 (
226 new DTRMParticle(mesh_, is, true)
227 );
228 }
229 };
230
231
232 // Access
233
234 //- Return const access to the initial position
235 inline const point& p0() const;
236
237 //- Return const access to the target position
238 inline const point& p1() const;
239
240 //- Return const access to the initial intensity
241 inline scalar I0() const;
242
243 //- Return const access to the current intensity
244 inline scalar I() const;
245
246 //- Return const access dA
247 inline scalar dA() const;
248
249
250 // Edit
251
252 //- Return access to the target position
253 inline point& p1();
254
255 //- Return access to the initial intensity
256 inline scalar& I0();
257
258 //- Return access to the current intensity
259 inline scalar& I();
260
261 //- Return access to dA
262 inline scalar& dA();
263
264 //- Return access to reflectedId
265 inline label& reflectedId();
266
267
268 // Tracking
269
270 //- Move
271 bool move(Cloud<DTRMParticle>& , trackingData&, const scalar);
272
273
274 // Member Operators
275
276 //- Overridable function to handle the particle hitting a processorPatch
278 (
280 trackingData& td
281 );
282
283 //- Overridable function to handle the particle hitting a wallPatch
284 void hitWallPatch
285 (
287 trackingData& td
288 );
290 bool hitPatch
291 (
293 trackingData& td
294 );
295
296
297 // I-O
298
299 //- Write individual parcel properties to stream
300 void writeProperties
301 (
302 Ostream& os,
303 const wordRes& filters,
304 const word& delim,
305 const bool namesOnly = false
306 ) const;
307
308
309 // Ostream Operator
311 friend Ostream& operator<<(Ostream& os, const DTRMParticle& p);
312};
313
314
315// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316
317} // End namespace Foam
318
319// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320
321#include "DTRMParticleI.H"
322
323// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324
325#endif
326
327// ************************************************************************* //
Base cloud calls templated on particle type.
Definition: Cloud.H:68
Factory class to read-construct particles used for.
Definition: DTRMParticle.H:211
iNew(const polyMesh &mesh)
Definition: DTRMParticle.H:216
autoPtr< DTRMParticle > operator()(Istream &is) const
Definition: DTRMParticle.H:221
Class used to pass tracking data to the trackToFace function.
Definition: DTRMParticle.H:95
const UPtrList< reflectionModel > & reflection() const
Definition: DTRMParticleI.H:98
const labelField & relfectedCells() const
Definition: DTRMParticleI.H:91
const interpolationCell< scalar > & aInterp() const
Definition: DTRMParticleI.H:58
const interpolationCellPoint< vector > & nHatInterp() const
Definition: DTRMParticleI.H:85
const interpolationCell< scalar > & EInterp() const
Definition: DTRMParticleI.H:72
const interpolationCell< scalar > & TInterp() const
Definition: DTRMParticleI.H:79
const interpolationCell< scalar > & eInterp() const
Definition: DTRMParticleI.H:65
Discrete Transfer Radiation Model (DTRM) particle.
Definition: DTRMParticle.H:65
scalar dA() const
Return const access dA.
DTRMParticle(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
const point & p0() const
Return const access to the initial position.
const point & p1() const
Return const access to the target position.
AddToPropertyList(particle, " p0"+" p1"+" I0"+" I"+" dA"+" transmissiveId";)
String representation of properties.
DTRMParticle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &position, const vector &targetPosition, const scalar I, const scalar dA, const label transmissiveId)
Construct from components.
scalar I() const
Return const access to the current intensity.
friend Ostream & operator<<(Ostream &os, const DTRMParticle &p)
DTRMParticle(const DTRMParticle &p)
Construct as copy.
label & reflectedId()
Return access to reflectedId.
void hitProcessorPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a processorPatch.
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: DTRMParticle.H:149
bool hitPatch(Cloud< DTRMParticle > &, trackingData &td)
scalar I0() const
Return const access to the initial intensity.
bool move(Cloud< DTRMParticle > &, trackingData &, const scalar)
Move.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
void hitWallPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
Uses the cell value for any location within the cell.
Base particle class.
Definition: particle.H:79
vector position() const
Return current particle position.
Definition: particleI.H:314
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition: particleI.H:143
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
Namespace for radiation modelling.
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.