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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
27 
28 #include "thermoCoupleProbes.H"
29 #include "mathematicalConstants.H"
30 #include "constants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
42  addToRunTimeSelectionTable(functionObject, thermoCoupleProbes, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::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 exist (resume old calculation)
72  // or of it is new.
73  dictionary probeDict;
74  if (getDict(typeName, probeDict))
75  {
76  probeDict.readEntry("Tc", Ttc_);
77  }
78  else
79  {
80  Ttc_ = probes::sample(thermo_.T());
81  }
82 
83  // Note: can only create the solver once all samples have been found
84  // - the number of samples is used to set the size of the ODE system
85  odeSolver_ = ODESolver::New(*this, dict);
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
96 {
97  return this->size();
98 }
99 
100 
102 (
103  const scalar x,
104  const scalarField& y,
105  scalarField& dydx
106 ) const
107 {
108  scalarField G(y.size(), Zero);
109  scalarField Tc(y.size(), Zero);
110  scalarField Uc(y.size(), Zero);
111  scalarField rhoc(y.size(), Zero);
112  scalarField muc(y.size(), Zero);
113  scalarField Cpc(y.size(), Zero);
114  scalarField kappac(y.size(), Zero);
115 
116  if (radiationFieldName_ != "none")
117  {
118  G = sample(mesh_.lookupObject<volScalarField>(radiationFieldName_));
119  }
120 
121  Tc = probes::sample(thermo_.T());
122 
123  Uc = mag(this->sample(mesh_.lookupObject<volVectorField>(UName_)));
124 
125  rhoc = this->sample(thermo_.rho()());
126  kappac = this->sample(thermo_.kappa()());
127  muc = this->sample(thermo_.mu()());
128  Cpc = this->sample(thermo_.Cp()());
129 
130  scalarField Re(rhoc*Uc*d_/muc);
131  scalarField Pr(Cpc*muc/kappac);
132  Pr = max(ROOTVSMALL, Pr);
133  scalarField Nu(2.0 + (0.4*sqrt(Re) + 0.06*pow(Re, 2.0/3.0))*pow(Pr, 0.4));
134  scalarField htc(Nu*kappac/d_);
135 
136  const scalar sigma = physicoChemical::sigma.value();
137 
138  scalar area = 4*constant::mathematical::pi*sqr(0.5*d_);
139  scalar volume = (4.0/3.0)*constant::mathematical::pi*pow3(0.5*d_);
140 
141  dydx =
142  (epsilon_*(G/4 - sigma*pow4(y))*area + htc*(Tc - y)*area)
143  / (rho_*Cp_*volume);
144 }
145 
146 
148 (
149  const scalar x,
150  const scalarField& y,
151  scalarField& dfdt,
152  scalarSquareMatrix& dfdy
153 ) const
154 {
155  derivatives(x, y, dfdt);
156 
157  const label n = nEqns();
158 
159  for (label i=0; i<n; i++)
160  {
161  for (label j=0; j<n; j++)
162  {
163  dfdy(i, j) = 0.0;
164  }
165  }
166 }
167 
168 
170 {
171  if (this->size())
172  {
173  sampleAndWrite<scalar>(thermo_.T());
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 (this->size())
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 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
runTime
engineTime & runTime
Definition: createEngineTime.H:13
mathematicalConstants.H
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::thermoCoupleProbes::derivatives
virtual void derivatives(const scalar x, const scalarField &y, scalarField &dydx) const
Calculate the derivatives in dydx.
Definition: thermoCoupleProbes.C:102
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::ODESolver::New
static autoPtr< ODESolver > New(const ODESystem &ode, const dictionary &dict)
Select null constructed.
Definition: ODESolverNew.C:34
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
Foam::functionObjects::thermoCoupleProbes::nEqns
virtual label nEqns() const
Number of ODE's to solve.
Definition: thermoCoupleProbes.C:95
Foam::fieldTypes::volume
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::functionObjects::thermoCoupleProbes::execute
virtual bool execute()
Execute, currently does nothing.
Definition: thermoCoupleProbes.C:185
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::probes::sample
tmp< Field< Type > > sample(const GeometricField< Type, fvPatchField, volMesh > &) const
Sample a volume field at all locations.
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::probes::read
virtual bool read(const dictionary &)
Read the probes.
Definition: probes.C:317
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::functionObjects::thermoCoupleProbes::~thermoCoupleProbes
virtual ~thermoCoupleProbes()
Destructor.
Definition: thermoCoupleProbes.C:91
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::functionObjects::thermoCoupleProbes::write
virtual bool write()
Sample and write.
Definition: thermoCoupleProbes.C:169
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::probes
Set of locations to sample.
Definition: probes.H:110
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::functionObjects::thermoCoupleProbes
Sample probe for temperature using a thermocouple.
Definition: thermoCoupleProbes.H:101
constants.H
Foam::SquareMatrix< scalar >
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::ODESystem
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:155
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:708
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::functionObjects::thermoCoupleProbes::jacobian
virtual void jacobian(const scalar t, const scalarField &y, scalarField &dfdt, scalarSquareMatrix &dfdy) const
Calculate the Jacobian of the system.
Definition: thermoCoupleProbes.C:148
Foam::functionObjects::thermoCoupleProbes::read
virtual bool read(const dictionary &)
Read.
Definition: thermoCoupleProbes.C:199
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
y
scalar y
Definition: LISASMDCalcMethod1.H:14
thermoCoupleProbes.H