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 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 
47  (
48  option,
49  solidificationMeltingSource,
50  dictionary
51  );
52  }
53 }
54 
55 const Foam::Enum
56 <
58 >
60 ({
61  { thermoMode::mdThermo, "thermo" },
62  { thermoMode::mdLookup, "lookup" },
63 });
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
69 Foam::fv::solidificationMeltingSource::Cp() const
70 {
71  switch (mode_)
72  {
73  case mdThermo:
74  {
75  const basicThermo& thermo =
77 
78  return thermo.Cp();
79  break;
80  }
81  case mdLookup:
82  {
83  if (CpName_ == "CpRef")
84  {
85  const scalar CpRef = coeffs_.get<scalar>("CpRef");
86 
88  (
89  IOobject
90  (
91  name_ + ":Cp",
92  mesh_.time().timeName(),
93  mesh_,
96  ),
97  mesh_,
99  (
100  "Cp",
102  CpRef
103  ),
104  extrapolatedCalculatedFvPatchScalarField::typeName
105  );
106  }
107  else
108  {
109  return mesh_.lookupObject<volScalarField>(CpName_);
110  }
111 
112  break;
113  }
114  default:
115  {
117  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
118  << abort(FatalError);
119  }
120  }
121 
122  return nullptr;
123 }
124 
125 
126 void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
127 {
128  if (curTimeIndex_ == mesh_.time().timeIndex())
129  {
130  return;
131  }
132 
133  if (debug)
134  {
135  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
136  }
137 
138  // update old time alpha1 field
139  alpha1_.oldTime();
140 
141  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
142 
143  forAll(cells_, i)
144  {
145  label celli = cells_[i];
146 
147  scalar Tc = T[celli];
148  scalar Cpc = Cp[celli];
149  scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
150 
151  alpha1_[celli] = max(0, min(alpha1New, 1));
152  deltaT_[i] = Tc - Tmelt_;
153  }
154 
155  alpha1_.correctBoundaryConditions();
156 
157  curTimeIndex_ = mesh_.time().timeIndex();
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162 
163 Foam::fv::solidificationMeltingSource::solidificationMeltingSource
164 (
165  const word& sourceName,
166  const word& modelType,
167  const dictionary& dict,
168  const fvMesh& mesh
169 )
170 :
171  cellSetOption(sourceName, modelType, dict, mesh),
172  Tmelt_(coeffs_.get<scalar>("Tmelt")),
173  L_(coeffs_.get<scalar>("L")),
174  relax_(coeffs_.lookupOrDefault("relax", 0.9)),
175  mode_(thermoModeTypeNames_.get("thermoMode", coeffs_)),
176  rhoRef_(coeffs_.get<scalar>("rhoRef")),
177  TName_(coeffs_.lookupOrDefault<word>("T", "T")),
178  CpName_(coeffs_.lookupOrDefault<word>("Cp", "Cp")),
179  UName_(coeffs_.lookupOrDefault<word>("U", "U")),
180  phiName_(coeffs_.lookupOrDefault<word>("phi", "phi")),
181  Cu_(coeffs_.lookupOrDefault<scalar>("Cu", 100000)),
182  q_(coeffs_.lookupOrDefault("q", 0.001)),
183  beta_(coeffs_.get<scalar>("beta")),
184  alpha1_
185  (
186  IOobject
187  (
188  name_ + ":alpha1",
189  mesh.time().timeName(),
190  mesh,
193  ),
194  mesh,
196  zeroGradientFvPatchScalarField::typeName
197  ),
198  curTimeIndex_(-1),
199  deltaT_(cells_.size(), 0)
200 {
201  fieldNames_.setSize(2);
202  fieldNames_[0] = UName_;
203 
204  switch (mode_)
205  {
206  case mdThermo:
207  {
208  const basicThermo& thermo =
209  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
210 
211  fieldNames_[1] = thermo.he().name();
212  break;
213  }
214  case mdLookup:
215  {
216  fieldNames_[1] = TName_;
217  break;
218  }
219  default:
220  {
222  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
223  << abort(FatalError);
224  }
225  }
226 
227  applied_.setSize(fieldNames_.size(), false);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232 
234 (
235  fvMatrix<scalar>& eqn,
236  const label fieldi
237 )
238 {
239  apply(geometricOneField(), eqn);
240 }
241 
242 
244 (
245  const volScalarField& rho,
246  fvMatrix<scalar>& eqn,
247  const label fieldi
248 )
249 {
250  apply(rho, eqn);
251 }
252 
253 
255 (
256  fvMatrix<vector>& eqn,
257  const label fieldi
258 )
259 {
260  if (debug)
261  {
262  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
263  }
264 
265  const volScalarField Cp(this->Cp());
266 
267  update(Cp);
268 
269  const vector& g = meshObjects::gravity::New(mesh_.time()).value();
270 
271  scalarField& Sp = eqn.diag();
272  vectorField& Su = eqn.source();
273  const scalarField& V = mesh_.V();
274 
275  forAll(cells_, i)
276  {
277  label celli = cells_[i];
278 
279  scalar Vc = V[celli];
280  scalar alpha1c = alpha1_[celli];
281 
282  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
283  vector Sb = rhoRef_*g*beta_*deltaT_[i];
284 
285  Sp[celli] += Vc*S;
286  Su[celli] += Vc*Sb;
287  }
288 }
289 
290 
292 (
293  const volScalarField& rho,
294  fvMatrix<vector>& eqn,
295  const label fieldi
296 )
297 {
298  // Momentum source uses a Boussinesq approximation - redirect
299  addSup(eqn, fieldi);
300 }
301 
302 
304 {
306  {
307  coeffs_.readEntry("Tmelt", Tmelt_);
308  coeffs_.readEntry("L", L_);
309 
310  coeffs_.readIfPresent("relax", relax_);
311 
312  thermoModeTypeNames_.readEntry("thermoMode", coeffs_, mode_);
313 
314  coeffs_.readEntry("rhoRef", rhoRef_);
315  coeffs_.readIfPresent("T", TName_);
316  coeffs_.readIfPresent("U", UName_);
317 
318  coeffs_.readIfPresent("Cu", Cu_);
319  coeffs_.readIfPresent("q", q_);
320 
321  coeffs_.readIfPresent("beta", beta_);
322 
323  return true;
324  }
325 
326  return false;
327 }
328 
329 
330 // ************************************************************************* //
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:51
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
Cell-set options abtract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:72
update
mesh update()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
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:54
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:337
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:76
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
rho
rho
Definition: readInitialConditions.H:96
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
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
Definition: solidificationMeltingSource.H:197
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:88
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:234
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_
Definition: solidificationMeltingSource.H:203
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:84
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:137
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
Foam::fv::solidificationMeltingSource::mdThermo
Definition: solidificationMeltingSource.H:199
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
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:303
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:200
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:295
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:246
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