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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::DTRMParticle
28 
29 Description
30  Discrete Transfer Radiation Model (DTRM) particle
31 
32 SourceFiles
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"
46 #include "interpolationCellPoint.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // Forward Declarations
54 class DTRMParticle;
55 
56 Ostream& operator<<(Ostream&, const DTRMParticle&);
57 
58 using namespace Foam::radiation;
59 
60 /*---------------------------------------------------------------------------*\
61  Class DTRMParticle Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 class 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 
89 public:
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_;
103  const interpolationCell<scalar>& TInterp_;
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  (
122  Cloud<DTRMParticle>& spc,
123  const interpolationCell<scalar>& aInterp,
124  const interpolationCell<scalar>& eInterp,
125  const interpolationCell<scalar>& EInterp,
126  const interpolationCell<scalar>& TInterp,
127  const interpolationCellPoint<vector>& nHatInterp,
128  const labelField&,
130  volScalarField& Q
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,
174  const vector& position,
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,
186  const barycentric& coordinates,
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
207  DTRMParticle(const DTRMParticle& p);
208 
209 
210  //- Factory class to read-construct particles used for
211  // parallel transfer
212  class iNew
213  {
214  const polyMesh& mesh_;
215 
216  public:
217 
218  iNew(const polyMesh& mesh)
219  :
220  mesh_(mesh)
221  {}
222 
223  autoPtr<DTRMParticle> operator()(Istream& is) const
224  {
225  return autoPtr<DTRMParticle>
226  (
227  new DTRMParticle(mesh_, is, true)
228  );
229  }
230  };
231 
232 
233  // Access
234 
235  //- Return const access to the initial position
236  const point& p0() const noexcept { return p0_; }
237 
238  //- Return const access to the target position
239  const point& p1() const noexcept { return p1_; }
240 
241  //- Return const access to the initial intensity
242  scalar I0() const noexcept { return I0_; }
244  //- Return const access to the current intensity
245  scalar I() const noexcept { return I_; }
246 
247  //- Return const access dA
248  scalar dA() const noexcept { return dA_; }
250 
251  // Edit
252 
253  //- Return access to the target position
254  point& p1() noexcept { return p1_; }
255 
256  //- Return access to the initial intensity
257  scalar& I0() noexcept { return I0_; }
258 
259  //- Return access to the current intensity
260  scalar& I() noexcept { return I_; }
261 
262  //- Return access to dA
263  scalar& dA() noexcept { return dA_; }
264 
265 
266  // Tracking
267 
268  //- Move
269  bool move(Cloud<DTRMParticle>& , trackingData&, const scalar);
270 
271 
272  // Member Operators
273 
274  //- Overridable function to handle the particle hitting a processorPatch
275  void hitProcessorPatch
276  (
278  trackingData& td
279  );
280 
281  //- Overridable function to handle the particle hitting a wallPatch
282  void hitWallPatch
283  (
285  trackingData& td
286  );
287 
288  bool hitPatch
289  (
291  trackingData& td
292  );
293 
294 
295  // I-O
296 
297  //- Write individual parcel properties to stream
298  void writeProperties
299  (
300  Ostream& os,
301  const wordRes& filters,
302  const word& delim,
303  const bool namesOnly = false
304  ) const;
305 
306 
307  // Ostream Operator
308 
309  friend Ostream& operator<<(Ostream& os, const DTRMParticle& p);
310 };
311 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 } // End namespace Foam
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 #include "DTRMParticleI.H"
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #endif
324 
325 // ************************************************************************* //
Class used to pass tracking data to the trackToFace function.
Definition: DTRMParticle.H:103
Forwards and collection of common volume field types.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Base particle class.
Definition: particle.H:69
#define AddToPropertyList(ParcelType, str)
Add to existing static &#39;propertyList&#39; for particle properties.
dynamicFvMesh & mesh
static const Identity< scalar > I
Definition: Identity.H:100
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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Factory class to read-construct particles used for.
Definition: DTRMParticle.H:243
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
Base cloud calls templated on particle type.
Definition: Cloud.H:51
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
PtrList< coordinateSystem > coordinates(solidRegions.size())
vector point
Point is a vector.
Definition: point.H:37
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Namespace for radiation modelling.
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
Namespace for OpenFOAM.
Discrete Transfer Radiation Model (DTRM) particle.
Definition: DTRMParticle.H:59