ThermoParcelTrackingDataI.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-2021 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 template<class ParcelType>
27 template<class TrackCloudType>
29 (
30  const TrackCloudType& cloud
31 )
32 :
33  ParcelType::trackingData(cloud),
34  Cp_(cloud.carrierThermo().Cp()),
35  kappa_(cloud.carrierThermo().kappa()),
36  pInterp_
37  (
38  interpolation<scalar>::New
39  (
40  cloud.solution().interpolationSchemes(),
41  cloud.p()
42  )
43  ),
44  TInterp_
45  (
46  interpolation<scalar>::New
47  (
48  cloud.solution().interpolationSchemes(),
49  cloud.T()
50  )
51  ),
52  CpInterp_
53  (
54  interpolation<scalar>::New
55  (
56  cloud.solution().interpolationSchemes(),
57  Cp_
58  )
59  ),
60  kappaInterp_
61  (
62  interpolation<scalar>::New
63  (
64  cloud.solution().interpolationSchemes(),
65  kappa_
66  )
67  ),
68  GInterp_(nullptr),
69  pc_(Zero),
70  Tc_(Zero),
71  Cpc_(Zero)
72 {
73  if (cloud.radiation())
74  {
75  GInterp_.reset
76  (
78  (
79  cloud.solution().interpolationSchemes(),
80  cloud.mesh().objectRegistry::template
81  lookupObject<volScalarField>("G")
82  ).ptr()
83  );
84  }
85 }
86 
87 
88 template<class ParcelType>
89 inline const Foam::volScalarField&
91 {
92  return Cp_;
93 }
94 
95 
96 template<class ParcelType>
97 inline const Foam::volScalarField&
99 {
100  return kappa_;
101 }
102 
103 
104 template<class ParcelType>
107 {
108  return pInterp_();
109 }
110 
111 
112 template<class ParcelType>
115 {
116  return TInterp_();
117 }
118 
119 
120 template<class ParcelType>
123 {
124  return CpInterp_();
125 }
126 
127 
128 template<class ParcelType>
131 {
132  return kappaInterp_();
133 }
134 
135 
136 template<class ParcelType>
139 {
140  if (!GInterp_.valid())
141  {
143  << "Radiation G interpolation object not set"
144  << abort(FatalError);
145  }
146 
147  return GInterp_();
148 }
149 
150 
151 template<class ParcelType>
153 {
154  return pc_;
155 }
156 
157 
158 template<class ParcelType>
160 {
161  return pc_;
162 }
163 
164 
165 template<class ParcelType>
167 {
168  return Tc_;
169 }
170 
171 
172 template<class ParcelType>
174 {
175  return Tc_;
176 }
177 
178 
179 template<class ParcelType>
181 {
182  return Cpc_;
183 }
184 
185 
186 template<class ParcelType>
188 {
189  return Cpc_;
190 }
191 
192 
193 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
Generic GeometricField class.
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
trackingData(const TrackCloudType &cloud)
Construct from components.
const volScalarField & kappa() const
Return access to the locally stored carrier kappa field.
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
scalar Cpc() const
Return the continuous phase specific heat capacity.
scalar pc() const
Return the continuous phase pressure.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & pInterp() const
Return const access to the interpolator for continuous phase.
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
const volScalarField & Cp() const
Return access to the locally stored carrier Cp field.
scalar Tc() const
Return the continuous phase temperature.
scalar Cp_
Specific heat capacity [J/kg/K].
Definition: ThermoParcel.H:271
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
Abstract base class for interpolation.
Definition: interpolation.H:55
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
error FatalError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p