Loading...
Searching...
No Matches
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-2024 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 Foam_DTRMParticle_H
38#define Foam_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
52
53// Forward Declarations
54class DTRMParticle;
55
57
58using namespace Foam::radiation;
60/*---------------------------------------------------------------------------*\
61 Class DTRMParticle Declaration
62\*---------------------------------------------------------------------------*/
63
64class DTRMParticle
65:
66 public particle
67{
68 // Private data
69
70 //- Initial position
71 point p0_;
72
73 //- Target position
74 point p1_;
75
76 //- Initial radiation intensity [W/m2]
77 scalar I0_;
78
79 //- Radiation intensity [W/m2]
80 scalar I_;
81
82 //- Area of radiation
83 scalar dA_;
84
85 //- Trasnmissive index
86 label transmissiveId_;
87
88
89public:
90
91 friend class Cloud<DTRMParticle>;
92
93 //- Class used to pass tracking data to the trackToFace function
94 class trackingData
95 :
97 {
98 // Interpolators for continuous phase fields
99
100 const interpolationCell<scalar>& aInterp_;
101 const interpolationCell<scalar>& eInterp_;
102 const interpolationCell<scalar>& EInterp_;
104 const interpolationCellPoint<vector>& nHatInterp_;
105
106 //- Reflected cells
107 const labelField& relfectedCells_;
108
109 //- Ptr to reflectiom model
110 UPtrList<reflectionModel> reflection_;
111
112 //- Heat source term
113 volScalarField& Q_;
114
115
116 public:
117
118 // Constructors
119
120 inline trackingData
121 (
128 const labelField&,
131 );
132
133
134 // Member Functions
135
136 inline const interpolationCell<scalar>& aInterp() const;
137 inline const interpolationCell<scalar>& eInterp() const;
138 inline const interpolationCell<scalar>& EInterp() const;
139 inline const interpolationCell<scalar>& TInterp() const;
140 inline const interpolationCellPoint<vector>& nHatInterp() const;
141 inline const labelField& relfectedCells() const;
142 inline const UPtrList<reflectionModel>& reflection() const;
143
144 inline scalar& Q(label celli);
145 };
146
147
148 // Static Data Members
149
150 //- Size in bytes of the fields
151 static const std::size_t sizeofFields_;
152
153
154 //- String representation of properties
156 (
157 particle,
158 " p0"
159 + " p1"
160 + " I0"
161 + " I"
162 + " dA"
163 + " transmissiveId";
164 );
165
166
167 // Constructors
169 //- Construct from components, with searching for tetFace and
170 // tetPt unless disabled by doCellFacePt = false.
172 (
173 const polyMesh& mesh,
175 const vector& targetPosition,
176 const scalar I,
177 const label cellI,
178 const scalar dA,
179 const label transmissiveId
180 );
181
182 //- Construct from components
184 (
185 const polyMesh& mesh,
187 const label celli,
188 const label tetFacei,
189 const label tetPti,
190 const vector& position,
191 const vector& targetPosition,
192 const scalar I,
193 const scalar dA,
194 const label transmissiveId
195 );
196
197 //- Construct from Istream
199 (
200 const polyMesh& mesh,
201 Istream& is,
202 bool readFields = true,
203 bool newFormat = true
204 );
205
206 //- Construct as copy
208
209 //- Return a clone
210 virtual autoPtr<particle> clone() const
211 {
212 return particle::Clone(*this);
213 }
214
215 //- Factory class to read-construct particles (for parallel transfer)
216 class iNew
217 {
218 const polyMesh& mesh_;
219
220 public:
221
222 iNew(const polyMesh& mesh)
223 :
224 mesh_(mesh)
225 {}
226
228 {
230 (
231 new DTRMParticle(mesh_, is, true)
232 );
233 }
234 };
236
237 // Access
238
239 //- Return const access to the initial position
240 const point& p0() const noexcept { return p0_; }
241
242 //- Return const access to the target position
243 const point& p1() const noexcept { return p1_; }
244
245 //- Return const access to the initial intensity
246 scalar I0() const noexcept { return I0_; }
247
248 //- Return const access to the current intensity
249 scalar I() const noexcept { return I_; }
250
251 //- Return const access dA
252 scalar dA() const noexcept { return dA_; }
253
255 // Edit
256
257 //- Return access to the target position
258 point& p1() noexcept { return p1_; }
260 //- Return access to the initial intensity
261 scalar& I0() noexcept { return I0_; }
262
263 //- Return access to the current intensity
264 scalar& I() noexcept { return I_; }
265
266 //- Return access to dA
267 scalar& dA() noexcept { return dA_; }
268
269
270 // Tracking
271
272 //- Move
273 bool move(Cloud<DTRMParticle>& , trackingData&, const scalar);
275
276 // Member Operators
277
278 //- Overridable function to handle the particle hitting a processorPatch
280 (
283 );
285 //- Overridable function to handle the particle hitting a wallPatch
286 void hitWallPatch
287 (
290 );
291
292 bool hitPatch
293 (
296 );
297
298
299 // I-O
300
301 //- Write individual parcel properties to stream
303 (
304 Ostream& os,
305 const wordRes& filters,
306 const word& delim,
307 const bool namesOnly = false
308 ) const;
309
310
311 // Ostream Operator
313 friend Ostream& operator<<(Ostream& os, const DTRMParticle& p);
314};
315
316
317// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318
319} // End namespace Foam
320
321// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322
323#include "DTRMParticleI.H"
324
325// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326
327#endif
328
329// ************************************************************************* //
Base cloud calls templated on particle type.
Definition Cloud.H:64
iNew(const polyMesh &mesh)
autoPtr< DTRMParticle > operator()(Istream &is) const
Class used to pass tracking data to the trackToFace function.
const UPtrList< reflectionModel > & reflection() const
trackingData(Cloud< DTRMParticle > &spc, const interpolationCell< scalar > &aInterp, const interpolationCell< scalar > &eInterp, const interpolationCell< scalar > &EInterp, const interpolationCell< scalar > &TInterp, const interpolationCellPoint< vector > &nHatInterp, const labelField &, const UPtrList< reflectionModel > &, volScalarField &Q)
const labelField & relfectedCells() const
const interpolationCell< scalar > & aInterp() const
const interpolationCellPoint< vector > & nHatInterp() const
const interpolationCell< scalar > & EInterp() const
const interpolationCell< scalar > & TInterp() const
const interpolationCell< scalar > & eInterp() const
Discrete Transfer Radiation Model (DTRM) particle.
DTRMParticle(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
scalar dA() const noexcept
Return const access dA.
point & p1() noexcept
Return access to the target position.
virtual autoPtr< particle > clone() const
Return a clone.
AddToPropertyList(particle, " p0"+" p1"+" I0"+" I"+" dA"+" transmissiveId";)
String representation of properties.
scalar I() const noexcept
Return const access to the current intensity.
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.
DTRMParticle(const polyMesh &mesh, const vector &position, const vector &targetPosition, const scalar I, const label cellI, const scalar dA, const label transmissiveId)
Construct from components, with searching for tetFace and.
const point & p0() const noexcept
Return const access to the initial position.
friend Ostream & operator<<(Ostream &os, const DTRMParticle &p)
scalar I0() const noexcept
Return const access to the initial intensity.
DTRMParticle(const DTRMParticle &p)
Construct as copy.
const point & p1() const noexcept
Return const access to the target position.
scalar & I0() noexcept
Return access to the initial intensity.
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.
bool hitPatch(Cloud< DTRMParticle > &, trackingData &td)
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.
scalar & I() noexcept
Return access to the current intensity.
scalar & dA() noexcept
Return access to dA.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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.
vector position() const
Return current particle position.
Definition particleI.H:283
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition particle.H:552
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for radiation modelling.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
Forwards and collection of common volume field types.