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-2020 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;
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace functionObjects
40{
43}
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
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 * * * * * * * * * * * * * * * //
89
91{
92 return this->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 if (radiationFieldName_ != "none")
112 {
113 G = sample(mesh_.lookupObject<volScalarField>(radiationFieldName_));
114 }
115
116 Tc = probes::sample(thermo_.T());
117
118 Uc = mag(this->sample(mesh_.lookupObject<volVectorField>(UName_)));
119
120 rhoc = this->sample(thermo_.rho()());
121 kappac = this->sample(thermo_.kappa()());
122 muc = this->sample(thermo_.mu()());
123 Cpc = this->sample(thermo_.Cp()());
124
125 scalarField Re(rhoc*Uc*d_/muc);
126 scalarField Pr(Cpc*muc/kappac);
127 Pr = max(ROOTVSMALL, Pr);
128 scalarField Nu(2.0 + (0.4*sqrt(Re) + 0.06*pow(Re, 2.0/3.0))*pow(Pr, 0.4));
129 scalarField htc(Nu*kappac/d_);
130
131 const scalar sigma = physicoChemical::sigma.value();
132
133 scalar area = 4*constant::mathematical::pi*sqr(0.5*d_);
134 scalar volume = (4.0/3.0)*constant::mathematical::pi*pow3(0.5*d_);
135
136 dydx =
137 (epsilon_*(G/4 - sigma*pow4(y))*area + htc*(Tc - y)*area)
138 / (rho_*Cp_*volume);
139}
140
141
143(
144 const scalar x,
145 const scalarField& y,
146 scalarField& dfdt,
148) const
149{
150 derivatives(x, y, dfdt);
151
152 const label n = nEqns();
153
154 for (label i=0; i<n; i++)
155 {
156 for (label j=0; j<n; j++)
157 {
158 dfdy(i, j) = 0.0;
159 }
160 }
161}
162
163
165{
166 if (!pointField::empty())
167 {
168 (void) prepare(ACTION_WRITE);
169
170 const auto& Tfield = thermo_.T();
171 writeValues(Tfield.name(), Ttc_, Tfield.time().timeOutputValue());
172
173 dictionary probeDict;
174 probeDict.add("Tc", Ttc_);
175 setProperty(typeName, probeDict);
176 return true;
177 }
178
179 return false;
180}
181
182
184{
185 if (!pointField::empty())
186 {
187 scalar dt = mesh_.time().deltaTValue();
188 scalar t = mesh_.time().value();
189 odeSolver_->solve(t, t+dt, Ttc_, dt);
190 return true;
191 }
192
193 return false;
194}
195
196
198{
199 if (probes::read(dict))
200 {
201 dict.readEntry("rho", rho_);
202 dict.readEntry("Cp", Cp_);
203 dict.readEntry("d", d_);
204 dict.readEntry("epsilon", epsilon_);
205 return true;
206 }
207
208 return false;
209}
210
211
212// ************************************************************************* //
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 const char *const typeName
Typename for Field.
Definition: FieldBase.H:59
Minimal example by using system/controlDict.functions:
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:50
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:622
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const Type & value() const
Return const reference to value.
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Abstract base-class for Time/database function objects.
Computes the magnitude of an input field.
Definition: mag.H:153
Computes the power of an input volScalarField.
Definition: pow.H:217
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
bool getDict(const word &entryName, dictionary &dict) const
Set dictionary, return true if set.
Sample probe for temperature using a thermocouple.
const fluidThermo & thermo_
Fluid thermo reference.
virtual void jacobian(const scalar t, const scalarField &y, scalarField &dfdt, scalarSquareMatrix &dfdy) const
Calculate the Jacobian of the system.
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. Evaluates the ODESolver.
scalarField Ttc_
Cached thermocouple temperature.
virtual bool write()
Sample and write.
virtual bool read(const dictionary &)
Read.
Set of locations to sample.
Definition: probes.H:166
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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.
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:47
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
dimensionedScalar Pr("Pr", dimless, laminarTransport)