Loading...
Searching...
No Matches
thermoCoupleProbes.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) 2016-2025 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
26\*---------------------------------------------------------------------------*/
27
28#include "thermoCoupleProbes.H"
30#include "constants.H"
32
33using namespace Foam::constant;
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace functionObjects
40{
43}
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49Foam::functionObjects::thermoCoupleProbes::thermoCoupleProbes
50(
51 const word& name,
52 const Time& runTime,
53 const dictionary& dict,
54 const bool loadFromFiles,
55 const bool readFields
56)
57:
58 probes(name, runTime, dict, loadFromFiles, false),
59 ODESystem(),
60 UName_(dict.getOrDefault<word>("U", "U")),
61 radiationFieldName_(dict.get<word>("radiationField")),
62 thermo_(mesh_.lookupObject<fluidThermo>(basicThermo::dictName)),
63 odeSolver_(nullptr),
64 Ttc_()
65{
66 if (readFields)
67 {
68 read(dict);
69 }
70
71 // Check if the property exists (resume old calculation) or is new
72 dictionary probeDict;
73 if (getDict(typeName, probeDict))
74 {
75 probeDict.readEntry("Tc", Ttc_);
76 }
77 else
78 {
80 }
81
82 // Note: can only create the solver once all samples have been found
83 // - the number of samples is used to set the size of the ODE system
85}
86
87
88// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91{
92 return probeModel().size();
93}
94
95
97(
98 const scalar x,
99 const scalarField& y,
100 scalarField& dydx
101) const
102{
103 scalarField G(y.size(), Zero);
104 scalarField Tc(y.size(), Zero);
105 scalarField Uc(y.size(), Zero);
106 scalarField rhoc(y.size(), Zero);
107 scalarField muc(y.size(), Zero);
108 scalarField Cpc(y.size(), Zero);
109 scalarField kappac(y.size(), Zero);
110
111 const auto& p = probeModel();
112
113 if (radiationFieldName_ != "none")
114 {
115 G = p.sample(mesh_.lookupObject<volScalarField>(radiationFieldName_));
116 }
117
118 Tc = p.sample(thermo_.T());
119
120 Uc = mag(p.sample(mesh_.lookupObject<volVectorField>(UName_)));
121
122 rhoc = p.sample(thermo_.rho()());
123 kappac = p.sample(thermo_.kappa()());
124 muc = p.sample(thermo_.mu()());
125 Cpc = p.sample(thermo_.Cp()());
126
127 scalarField Re(rhoc*Uc*d_/muc);
128 scalarField Pr(Cpc*muc/kappac);
129 Pr = max(ROOTVSMALL, Pr);
130 scalarField Nu(2.0 + (0.4*sqrt(Re) + 0.06*pow(Re, 2.0/3.0))*pow(Pr, 0.4));
131 scalarField htc(Nu*kappac/d_);
132
133 const scalar sigma = physicoChemical::sigma.value();
134
135 scalar area = 4*constant::mathematical::pi*sqr(0.5*d_);
136 scalar volume = (4.0/3.0)*constant::mathematical::pi*pow3(0.5*d_);
138 dydx =
139 (epsilon_*(G/4 - sigma*pow4(y))*area + htc*(Tc - y)*area)
140 / (rho_*Cp_*volume);
141}
142
143
145(
146 const scalar x,
147 const scalarField& y,
148 scalarField& dfdt,
150) const
151{
152 derivatives(x, y, dfdt);
153
154 const label n = nEqns();
155
156 for (label i=0; i<n; i++)
157 {
158 for (label j=0; j<n; j++)
160 dfdy(i, j) = 0.0;
161 }
162 }
163}
164
165
167{
168 if (!probeModel().empty())
169 {
170 (void) prepare(ACTION_WRITE);
171
172 const auto& Tfield = thermo_.T();
173 writeValues(Tfield.name(), Ttc_, Tfield.time().timeOutputValue());
174
175 dictionary probeDict;
176 probeDict.add("Tc", Ttc_);
177 setProperty(typeName, probeDict);
178 return true;
179 }
180
181 return false;
182}
183
184
186{
187 if (!probeModel().empty())
188 {
189 scalar dt = mesh_.time().deltaTValue();
190 scalar t = mesh_.time().value();
191 odeSolver_->solve(t, t+dt, Ttc_, dt);
192 return true;
193 }
194
195 return false;
196}
197
198
200{
201 if (probes::read(dict))
202 {
203 dict.readEntry("rho", rho_);
204 dict.readEntry("Cp", Cp_);
205 dict.readEntry("d", d_);
206 dict.readEntry("epsilon", epsilon_);
207 return true;
208 }
209
210 return false;
211}
212
213
214// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static autoPtr< ODESolver > New(const ODESystem &ode, const dictionary &dict)
Select null constructed.
Abstract base class for the systems of ordinary differential equations.
Definition ODESystem.H:44
ODESystem()
Construct null.
Definition ODESystem.H:53
label prepare(unsigned request)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
Abstract base-class for Time/database function objects.
const fvMesh & mesh_
Reference to the fvMesh.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool getDict(const word &entryName, dictionary &dict) const
Set dictionary, return true if set.
void setProperty(const word &entryName, const Type &value)
Add generic property.
Sample probe for temperature using a thermocouple.
const fluidThermo & thermo_
Fluid thermo reference.
scalar Cp_
Thermocouple heat capacity.
word UName_
Name of the velocity field.
virtual void jacobian(const scalar t, const scalarField &y, scalarField &dfdt, scalarSquareMatrix &dfdy) const
Calculate the Jacobian of the system.
scalar epsilon_
Thermocouple emissivity.
word radiationFieldName_
Name of the incident radiation field.
autoPtr< ODESolver > odeSolver_
ODESolver.
virtual void derivatives(const scalar x, const scalarField &y, scalarField &dydx) const
Calculate the derivatives in dydx.
virtual label nEqns() const
Number of ODE's to solve.
virtual bool execute()
Execute the function-object operations.
scalarField Ttc_
Cached thermocouple temperature.
virtual bool write()
Write the function-object results.
virtual bool read(const dictionary &)
Read the function-object dictionary.
tmp< Field< Type > > sample(const VolumeField< Type > &) const
Sample a volume field at all locations.
Base class for sampling fields at specified internal and boundary locations.
Definition probeModel.H:49
Function object to sample fields at specified internal-mesh locations and write to file.
Definition probes.H:169
const internalFieldProbe & probeModel() const noexcept
Bring Base::probeModel into this class's public scope.
Definition Probes.H:250
virtual bool read(const dictionary &)
Read the settings from the dictionary.
Definition probes.C:67
probes(const word &name, const Time &runTime, const dictionary &dict, const bool loadFromFiles=false, const bool readFields=true)
Construct from Time and dictionary.
Definition probes.C:41
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
engineTime & runTime
const word dictName("faMeshDefinition")
constexpr scalar pi(M_PI)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
SquareMatrix< scalar > scalarSquareMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
dimensionedScalar Pr("Pr", dimless, laminarTransport)