solidificationMeltingSource.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) 2014-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvMatrices.H"
31 #include "basicThermo.H"
32 #include "gravityMeshObject.H"
36 #include "geometricOneField.H"
37 
38 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace fv
43 {
44  defineTypeNameAndDebug(solidificationMeltingSource, 0);
45  addToRunTimeSelectionTable(option, solidificationMeltingSource, dictionary);
46 }
47 }
48 
49 const Foam::Enum
50 <
52 >
54 ({
55  { thermoMode::mdThermo, "thermo" },
56  { thermoMode::mdLookup, "lookup" },
57 });
58 
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
63 Foam::fv::solidificationMeltingSource::Cp() const
64 {
65  switch (mode_)
66  {
67  case mdThermo:
68  {
69  const auto& thermo =
71 
72  return thermo.Cp();
73  break;
74  }
75  case mdLookup:
76  {
77  if (CpName_ == "CpRef")
78  {
79  const scalar CpRef = coeffs_.get<scalar>("CpRef");
80 
82  (
83  IOobject
84  (
85  name_ + ":Cp",
86  mesh_.time().timeName(),
87  mesh_,
90  ),
91  mesh_,
93  (
94  "Cp",
96  CpRef
97  ),
98  extrapolatedCalculatedFvPatchScalarField::typeName
99  );
100  }
101  else
102  {
103  return mesh_.lookupObject<volScalarField>(CpName_);
104  }
105 
106  break;
107  }
108  default:
109  {
111  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
112  << abort(FatalError);
113  }
114  }
115 
116  return nullptr;
117 }
118 
119 
120 void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
121 {
122  if (curTimeIndex_ == mesh_.time().timeIndex())
123  {
124  return;
125  }
126 
127  if (debug)
128  {
129  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
130  }
131 
132  // update old time alpha1 field
133  alpha1_.oldTime();
134 
135  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
136 
137  forAll(cells_, i)
138  {
139  label celli = cells_[i];
140 
141  scalar Tc = T[celli];
142  scalar Cpc = Cp[celli];
143  scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
144 
145  alpha1_[celli] = max(0, min(alpha1New, 1));
146  deltaT_[i] = Tc - Tmelt_;
147  }
148 
149  alpha1_.correctBoundaryConditions();
150 
151  curTimeIndex_ = mesh_.time().timeIndex();
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
156 
158 (
159  const word& sourceName,
160  const word& modelType,
161  const dictionary& dict,
162  const fvMesh& mesh
163 )
164 :
165  cellSetOption(sourceName, modelType, dict, mesh),
166  Tmelt_(coeffs_.get<scalar>("Tmelt")),
167  L_(coeffs_.get<scalar>("L")),
168  relax_(coeffs_.getOrDefault<scalar>("relax", 0.9)),
169  mode_(thermoModeTypeNames_.get("thermoMode", coeffs_)),
170  rhoRef_(coeffs_.get<scalar>("rhoRef")),
171  TName_(coeffs_.getOrDefault<word>("T", "T")),
172  CpName_(coeffs_.getOrDefault<word>("Cp", "Cp")),
173  UName_(coeffs_.getOrDefault<word>("U", "U")),
174  phiName_(coeffs_.getOrDefault<word>("phi", "phi")),
175  Cu_(coeffs_.getOrDefault<scalar>("Cu", 100000)),
176  q_(coeffs_.getOrDefault<scalar>("q", 0.001)),
177  beta_(coeffs_.get<scalar>("beta")),
178  alpha1_
179  (
180  IOobject
181  (
182  name_ + ":alpha1",
183  mesh.time().timeName(),
184  mesh,
187  ),
188  mesh,
190  zeroGradientFvPatchScalarField::typeName
191  ),
192  curTimeIndex_(-1),
193  deltaT_(cells_.size(), 0)
194 {
195  fieldNames_.setSize(2);
196  fieldNames_[0] = UName_;
197 
198  switch (mode_)
199  {
200  case mdThermo:
201  {
202  const auto& thermo =
203  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
204 
205  fieldNames_[1] = thermo.he().name();
206  break;
207  }
208  case mdLookup:
209  {
210  fieldNames_[1] = TName_;
211  break;
212  }
213  default:
214  {
216  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
217  << abort(FatalError);
218  }
219  }
220 
221  applied_.setSize(fieldNames_.size(), false);
222 }
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
228 (
229  fvMatrix<scalar>& eqn,
230  const label fieldi
231 )
232 {
233  apply(geometricOneField(), eqn);
234 }
235 
236 
238 (
239  const volScalarField& rho,
240  fvMatrix<scalar>& eqn,
241  const label fieldi
242 )
243 {
244  apply(rho, eqn);
245 }
246 
247 
249 (
250  fvMatrix<vector>& eqn,
251  const label fieldi
252 )
253 {
254  if (debug)
255  {
256  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
257  }
258 
259  const volScalarField Cp(this->Cp());
260 
261  update(Cp);
262 
263  const vector& g = meshObjects::gravity::New(mesh_.time()).value();
264 
265  scalarField& Sp = eqn.diag();
266  vectorField& Su = eqn.source();
267  const scalarField& V = mesh_.V();
268 
269  forAll(cells_, i)
270  {
271  label celli = cells_[i];
272 
273  scalar Vc = V[celli];
274  scalar alpha1c = alpha1_[celli];
275 
276  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
277  vector Sb(rhoRef_*g*beta_*deltaT_[i]);
278 
279  Sp[celli] += Vc*S;
280  Su[celli] += Vc*Sb;
281  }
282 }
283 
284 
286 (
287  const volScalarField& rho,
288  fvMatrix<vector>& eqn,
289  const label fieldi
290 )
291 {
292  // Momentum source uses a Boussinesq approximation - redirect
293  addSup(eqn, fieldi);
294 }
295 
296 
298 {
300  {
301  coeffs_.readEntry("Tmelt", Tmelt_);
302  coeffs_.readEntry("L", L_);
303 
304  coeffs_.readIfPresent("relax", relax_);
305 
306  thermoModeTypeNames_.readEntry("thermoMode", coeffs_, mode_);
307 
308  coeffs_.readEntry("rhoRef", rhoRef_);
309  coeffs_.readIfPresent("T", TName_);
310  coeffs_.readIfPresent("U", UName_);
311 
312  coeffs_.readIfPresent("Cu", Cu_);
313  coeffs_.readIfPresent("q", q_);
314 
315  coeffs_.readEntry("beta", beta_);
316 
317  return true;
318  }
319 
320  return false;
321 }
322 
323 
324 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
solidificationMeltingSource.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
basicThermo.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Definition: cellSetOption.H:163
update
mesh update()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::fv::solidificationMeltingSource::solidificationMeltingSource
solidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: solidificationMeltingSource.C:158
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Sp
zeroField Sp
Definition: alphaSuSp.H:2
gravityMeshObject.H
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:133
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
rho
rho
Definition: readInitialConditions.H:88
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:287
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::fv::solidificationMeltingSource::thermoMode
thermoMode
Options for the thermo mode specification.
Definition: solidificationMeltingSource.H:260
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fv::option::coeffs_
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:145
Foam::fv::cellSetOption::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: cellSetOption.C:255
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::fv::solidificationMeltingSource::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
Definition: solidificationMeltingSource.C:228
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
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::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::fv::solidificationMeltingSource::thermoModeTypeNames_
static const Enum< thermoMode > thermoModeTypeNames_
Names for thermoMode.
Definition: solidificationMeltingSource.H:267
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
fv
labelList fv(nPoints)
Foam::fv::solidificationMeltingSource::mdThermo
Definition: solidificationMeltingSource.H:262
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
geometricOneField.H
Foam::fv::solidificationMeltingSource::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: solidificationMeltingSource.C:297
Foam::Vector< scalar >
Foam::apply
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Definition: parcelSelectionDetail.C:82
Foam::fv::solidificationMeltingSource::mdLookup
Definition: solidificationMeltingSource.H:263
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:297
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
zeroGradientFvPatchFields.H
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
extrapolatedCalculatedFvPatchFields.H